{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "# LETID - Outdoor Environments \n", "\n", "This is an example on how to model LETID progression in outdoor environments\n", "\n", "We can use the equations in this library to model LETID progression in a simulated outdoor environment, given that we have weather and system data. This example makes use of tools from the fabulous [pvlib](https://pvlib-python.readthedocs.io/en/stable/) library to calculate system irradiance and temperature, which we use to calculate progression in LETID states.\n", "\n", "This will illustrate the potential of \"Temporary Recovery\", i.e., the backwards transition of the LETID defect B->A that can take place with carrier injection at lower temperatures.\n", "\n", "\n", "**Requirements:**\n", "- `pvlib`, `pandas`, `numpy`, `matplotlib`\n", "\n", "**Objectives:**\n", "1. Use `pvlib` and provided weather files to set up a temperature and injection timeseries\n", "2. Define necessary solar cell device parameters\n", "3. Define necessary degradation parameters: degraded lifetime and defect states\n", "4. Run through timeseries, calculating defect states\n", "5. Calculate device degradation and plot\n" ] }, { "cell_type": "code", "execution_count": 1, "metadata": {}, "outputs": [], "source": [ "# if running on google colab, uncomment the next line and execute this cell to install the dependencies and prevent \"ModuleNotFoundError\" in later cells:\n", "# !pip install pvdeg==0.3.3" ] }, { "cell_type": "code", "execution_count": 2, "metadata": {}, "outputs": [], "source": [ "from pvdeg import letid, collection, utilities, DATA_DIR\n", "\n", "import pvlib\n", "import os\n", "import pandas as pd\n", "import numpy as np\n", "import matplotlib.pyplot as plt\n", "import pvdeg" ] }, { "cell_type": "code", "execution_count": 3, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Working on a Windows 10\n", "Python version 3.11.7 | packaged by Anaconda, Inc. | (main, Dec 15 2023, 18:05:47) [MSC v.1916 64 bit (AMD64)]\n", "Pandas version 2.2.0\n", "pvlib version 0.10.3\n", "pvdeg version 0.2.4.dev83+ge2ceab9.d20240422\n" ] } ], "source": [ "# This information helps with debugging and getting support :)\n", "import sys, platform\n", "print(\"Working on a \", platform.system(), platform.release())\n", "print(\"Python version \", sys.version)\n", "print(\"Pandas version \", pd.__version__)\n", "print(\"pvlib version \", pvlib.__version__)\n", "print(\"pvdeg version \", pvdeg.__version__)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "First, we'll use pvlib to create and run a model system, and use the irradiance, temperature, and operating point of that model to set up our LETID model\n", "For this example, we'll model a fixed latitude tilt system at NREL, in Golden, CO, USA, using [NSRDB](https://nsrdb.nrel.gov/) hourly PSM weather data, SAPM temperature models, and module and inverter models from the CEC database." ] }, { "cell_type": "code", "execution_count": 4, "metadata": {}, "outputs": [], "source": [ "# load weather and location data, use pvlib read_psm3 function with map_variables = True\n", "\n", "sam_file = 'psm3.csv'\n", "weather, meta = pvdeg.weather.read(os.path.join(DATA_DIR, sam_file), file_type='PSM3', map_variables = True)" ] }, { "cell_type": "code", "execution_count": 5, "metadata": {}, "outputs": [ { "data": { "text/html": [ "
\n", "\n", "\n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", "
YearMonthDayHourMinutednidhighitemp_airdew_pointwind_speedrelative_humiditypoa_globaltemp_celltemp_module
1999-01-01 00:30:00-07:001999110300.00.00.00.0-5.01.879.390.00.00.0
1999-01-01 01:30:00-07:001999111300.00.00.00.0-4.01.780.840.00.00.0
1999-01-01 02:30:00-07:001999112300.00.00.00.0-4.01.582.980.00.00.0
1999-01-01 03:30:00-07:001999113300.00.00.00.0-4.01.385.010.00.00.0
1999-01-01 04:30:00-07:001999114300.00.00.00.0-4.01.385.810.00.00.0
................................................
1999-12-31 19:30:00-07:001999123119300.00.00.00.0-3.00.983.630.00.00.0
1999-12-31 20:30:00-07:001999123120300.00.00.00.0-3.01.286.820.00.00.0
1999-12-31 21:30:00-07:001999123121300.00.00.00.0-4.01.683.780.00.00.0
1999-12-31 22:30:00-07:001999123122300.00.00.00.0-4.01.781.220.00.00.0
1999-12-31 23:30:00-07:001999123123300.00.00.00.0-5.01.879.430.00.00.0
\n", "

8760 rows × 15 columns

\n", "
" ], "text/plain": [ " Year Month Day Hour Minute dni dhi ghi \\\n", "1999-01-01 00:30:00-07:00 1999 1 1 0 30 0.0 0.0 0.0 \n", "1999-01-01 01:30:00-07:00 1999 1 1 1 30 0.0 0.0 0.0 \n", "1999-01-01 02:30:00-07:00 1999 1 1 2 30 0.0 0.0 0.0 \n", "1999-01-01 03:30:00-07:00 1999 1 1 3 30 0.0 0.0 0.0 \n", "1999-01-01 04:30:00-07:00 1999 1 1 4 30 0.0 0.0 0.0 \n", "... ... ... ... ... ... ... ... ... \n", "1999-12-31 19:30:00-07:00 1999 12 31 19 30 0.0 0.0 0.0 \n", "1999-12-31 20:30:00-07:00 1999 12 31 20 30 0.0 0.0 0.0 \n", "1999-12-31 21:30:00-07:00 1999 12 31 21 30 0.0 0.0 0.0 \n", "1999-12-31 22:30:00-07:00 1999 12 31 22 30 0.0 0.0 0.0 \n", "1999-12-31 23:30:00-07:00 1999 12 31 23 30 0.0 0.0 0.0 \n", "\n", " temp_air dew_point wind_speed relative_humidity \\\n", "1999-01-01 00:30:00-07:00 0.0 -5.0 1.8 79.39 \n", "1999-01-01 01:30:00-07:00 0.0 -4.0 1.7 80.84 \n", "1999-01-01 02:30:00-07:00 0.0 -4.0 1.5 82.98 \n", "1999-01-01 03:30:00-07:00 0.0 -4.0 1.3 85.01 \n", "1999-01-01 04:30:00-07:00 0.0 -4.0 1.3 85.81 \n", "... ... ... ... ... \n", "1999-12-31 19:30:00-07:00 0.0 -3.0 0.9 83.63 \n", "1999-12-31 20:30:00-07:00 0.0 -3.0 1.2 86.82 \n", "1999-12-31 21:30:00-07:00 0.0 -4.0 1.6 83.78 \n", "1999-12-31 22:30:00-07:00 0.0 -4.0 1.7 81.22 \n", "1999-12-31 23:30:00-07:00 0.0 -5.0 1.8 79.43 \n", "\n", " poa_global temp_cell temp_module \n", "1999-01-01 00:30:00-07:00 0.0 0.0 0.0 \n", "1999-01-01 01:30:00-07:00 0.0 0.0 0.0 \n", "1999-01-01 02:30:00-07:00 0.0 0.0 0.0 \n", "1999-01-01 03:30:00-07:00 0.0 0.0 0.0 \n", "1999-01-01 04:30:00-07:00 0.0 0.0 0.0 \n", "... ... ... ... \n", "1999-12-31 19:30:00-07:00 0.0 0.0 0.0 \n", "1999-12-31 20:30:00-07:00 0.0 0.0 0.0 \n", "1999-12-31 21:30:00-07:00 0.0 0.0 0.0 \n", "1999-12-31 22:30:00-07:00 0.0 0.0 0.0 \n", "1999-12-31 23:30:00-07:00 0.0 0.0 0.0 \n", "\n", "[8760 rows x 15 columns]" ] }, "execution_count": 5, "metadata": {}, "output_type": "execute_result" } ], "source": [ "weather" ] }, { "cell_type": "code", "execution_count": 6, "metadata": {}, "outputs": [], "source": [ "#if our weather file doesn't have precipitable water, calculate it with pvlib\n", "if not 'precipitable_water' in weather.columns:\n", " weather['precipitable_water'] = pvlib.atmosphere.gueymard94_pw(weather['temp_air'], weather['relative_humidity'])\n" ] }, { "cell_type": "code", "execution_count": 7, "metadata": {}, "outputs": [], "source": [ "# rename some columns for pvlib if they haven't been already\n", "weather.rename(columns = {'GHI':'ghi', 'DNI':'dni', 'DHI':'dhi', 'Temperature':'temp_air', 'Wind Speed':'wind_speed', 'Relative Humidity':'relative_humidity', 'Precipitable Water':'precipitable_water'}, inplace = True)\n", "weather = weather[['ghi', 'dni', 'dhi', 'temp_air', 'wind_speed', 'relative_humidity', 'precipitable_water']]" ] }, { "cell_type": "code", "execution_count": 8, "metadata": {}, "outputs": [ { "data": { "text/html": [ "
\n", "\n", "\n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", "
ghidnidhitemp_airwind_speedrelative_humidityprecipitable_water
1999-01-01 00:30:00-07:000.00.00.00.01.879.390.891869
1999-01-01 01:30:00-07:000.00.00.00.01.780.840.908158
1999-01-01 02:30:00-07:000.00.00.00.01.582.980.932199
1999-01-01 03:30:00-07:000.00.00.00.01.385.010.955004
1999-01-01 04:30:00-07:000.00.00.00.01.385.810.963991
........................
1999-12-31 19:30:00-07:000.00.00.00.00.983.630.939501
1999-12-31 20:30:00-07:000.00.00.00.01.286.820.975338
1999-12-31 21:30:00-07:000.00.00.00.01.683.780.941186
1999-12-31 22:30:00-07:000.00.00.00.01.781.220.912427
1999-12-31 23:30:00-07:000.00.00.00.01.879.430.892318
\n", "

8760 rows × 7 columns

\n", "
" ], "text/plain": [ " ghi dni dhi temp_air wind_speed \\\n", "1999-01-01 00:30:00-07:00 0.0 0.0 0.0 0.0 1.8 \n", "1999-01-01 01:30:00-07:00 0.0 0.0 0.0 0.0 1.7 \n", "1999-01-01 02:30:00-07:00 0.0 0.0 0.0 0.0 1.5 \n", "1999-01-01 03:30:00-07:00 0.0 0.0 0.0 0.0 1.3 \n", "1999-01-01 04:30:00-07:00 0.0 0.0 0.0 0.0 1.3 \n", "... ... ... ... ... ... \n", "1999-12-31 19:30:00-07:00 0.0 0.0 0.0 0.0 0.9 \n", "1999-12-31 20:30:00-07:00 0.0 0.0 0.0 0.0 1.2 \n", "1999-12-31 21:30:00-07:00 0.0 0.0 0.0 0.0 1.6 \n", "1999-12-31 22:30:00-07:00 0.0 0.0 0.0 0.0 1.7 \n", "1999-12-31 23:30:00-07:00 0.0 0.0 0.0 0.0 1.8 \n", "\n", " relative_humidity precipitable_water \n", "1999-01-01 00:30:00-07:00 79.39 0.891869 \n", "1999-01-01 01:30:00-07:00 80.84 0.908158 \n", "1999-01-01 02:30:00-07:00 82.98 0.932199 \n", "1999-01-01 03:30:00-07:00 85.01 0.955004 \n", "1999-01-01 04:30:00-07:00 85.81 0.963991 \n", "... ... ... \n", "1999-12-31 19:30:00-07:00 83.63 0.939501 \n", "1999-12-31 20:30:00-07:00 86.82 0.975338 \n", "1999-12-31 21:30:00-07:00 83.78 0.941186 \n", "1999-12-31 22:30:00-07:00 81.22 0.912427 \n", "1999-12-31 23:30:00-07:00 79.43 0.892318 \n", "\n", "[8760 rows x 7 columns]" ] }, "execution_count": 8, "metadata": {}, "output_type": "execute_result" } ], "source": [ "weather\n" ] }, { "cell_type": "code", "execution_count": 9, "metadata": {}, "outputs": [], "source": [ "# import pvlib stuff and pick a module and inverter. Choice of these things will slightly affect the pvlib results which we later use to calculate injection.\n", "# we'll use the SAPM temperature model open-rack glass/polymer coeffecients.\n", "\n", "from pvlib.temperature import TEMPERATURE_MODEL_PARAMETERS\n", "from pvlib.location import Location\n", "from pvlib.pvsystem import PVSystem\n", "from pvlib.modelchain import ModelChain\n", "\n", "cec_modules = pvlib.pvsystem.retrieve_sam('CECMod')\n", "cec_inverters = pvlib.pvsystem.retrieve_sam('cecinverter')\n", "\n", "cec_module = cec_modules['Jinko_Solar_Co___Ltd_JKM260P_60']\n", "cec_inverter = cec_inverters['ABB__ULTRA_750_TL_OUTD_1_US_690_x_y_z__690V_']\n", "\n", "temperature_model_parameters = TEMPERATURE_MODEL_PARAMETERS['sapm']['open_rack_glass_polymer']\n" ] }, { "cell_type": "code", "execution_count": 10, "metadata": {}, "outputs": [], "source": [ "# set up system in pvlib\n", "lat = meta['latitude']\n", "lon = meta['longitude']\n", "tz = meta['tz']\n", "elevation = meta['altitude']\n", "surface_tilt = lat # fixed, latitude tilt\n", "surface_azimuth = 180 # south-facing\n", "\n", "location = Location(lat, lon, tz, elevation, 'Golden, CO, USA')\n", "\n", "system = PVSystem(surface_tilt = surface_tilt, surface_azimuth = surface_azimuth,\n", " module_parameters = cec_module,\n", " inverter_parameters = cec_inverter,\n", " temperature_model_parameters = temperature_model_parameters,\n", " )" ] }, { "cell_type": "code", "execution_count": 11, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "ModelChain: \n", " name: None\n", " clearsky_model: ineichen\n", " transposition_model: haydavies\n", " solar_position_method: nrel_numpy\n", " airmass_model: kastenyoung1989\n", " dc_model: cec\n", " ac_model: sandia_inverter\n", " aoi_model: physical_aoi_loss\n", " spectral_model: first_solar_spectral_loss\n", " temperature_model: sapm_temp\n", " losses_model: no_extra_losses" ] }, "execution_count": 11, "metadata": {}, "output_type": "execute_result" } ], "source": [ "# create and run pvlib modelchain\n", "mc = ModelChain(system, location, aoi_model = \"physical\")\n", "mc.run_model(weather)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Set up timeseries\n", "In this example, injection is a function of both the operating point of the module (which we will assume is maximum power point) and irradiance. Maximum power point injection is equivalent to $(I_{sc}-I_{mp})/I_{sc}\\times Ee$, where $Ee$ is effective irradiance, the irradiance absorbed by the module's cells. We normalize it to 1-sun irradiance, 1000 $W/m^2$.\n", "\n", "We will use the irradiance, DC operating point, and cell temperature from the pvlib modelchain results." ] }, { "cell_type": "code", "execution_count": 12, "metadata": {}, "outputs": [], "source": [ "ee = mc.results.effective_irradiance\n", "#injection = (mc.results.dc['i_sc']-mc.results.dc['i_mp'])/(mc.results.dc['i_sc'])*(ee/1000)\n", "injection = letid.calc_injection_outdoors(mc.results)\n", "temperature = mc.results.cell_temperature\n", "\n", "timesteps = pd.DataFrame({'Temperature':temperature, 'Injection':injection}) # create a DataFrame with cell temperature and injection\n", "timesteps.reset_index(inplace = True) # reset the index so datetime is a column. I prefer integer indexing.\n", "timesteps.rename(columns = {'index':'Datetime'}, inplace = True)\n" ] }, { "cell_type": "code", "execution_count": 13, "metadata": {}, "outputs": [], "source": [ "# filter out times when injection is NaN, these won't progress LETID, and it'll make the calculations below run faster\n", "timesteps = timesteps[timesteps['Injection'].notnull()]\n", "timesteps.reset_index(inplace = True, drop = True)" ] }, { "cell_type": "code", "execution_count": 14, "metadata": {}, "outputs": [ { "data": { "text/html": [ "
\n", "\n", "\n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", "
DatetimeTemperatureInjection
01999-01-01 08:30:00-07:003.3487740.006239
11999-01-01 09:30:00-07:0012.3784390.031607
21999-01-01 10:30:00-07:0016.5950810.044220
31999-01-01 11:30:00-07:0017.4572400.043201
41999-01-01 12:30:00-07:006.3927830.007124
............
42961999-12-31 12:30:00-07:009.5155770.001151
42971999-12-31 13:30:00-07:0028.1197960.044499
42981999-12-31 14:30:00-07:0023.3146720.035130
42991999-12-31 15:30:00-07:0017.8905280.027891
43001999-12-31 16:30:00-07:004.5523650.001188
\n", "

4301 rows × 3 columns

\n", "
" ], "text/plain": [ " Datetime Temperature Injection\n", "0 1999-01-01 08:30:00-07:00 3.348774 0.006239\n", "1 1999-01-01 09:30:00-07:00 12.378439 0.031607\n", "2 1999-01-01 10:30:00-07:00 16.595081 0.044220\n", "3 1999-01-01 11:30:00-07:00 17.457240 0.043201\n", "4 1999-01-01 12:30:00-07:00 6.392783 0.007124\n", "... ... ... ...\n", "4296 1999-12-31 12:30:00-07:00 9.515577 0.001151\n", "4297 1999-12-31 13:30:00-07:00 28.119796 0.044499\n", "4298 1999-12-31 14:30:00-07:00 23.314672 0.035130\n", "4299 1999-12-31 15:30:00-07:00 17.890528 0.027891\n", "4300 1999-12-31 16:30:00-07:00 4.552365 0.001188\n", "\n", "[4301 rows x 3 columns]" ] }, "execution_count": 14, "metadata": {}, "output_type": "execute_result" } ], "source": [ "timesteps" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Device parameters\n", "To define a device, we need to define several important quantities about the device: wafer thickness (in $\\mu m$), rear surface recombination velocity (in cm/s), and cell area (in cm2)." ] }, { "cell_type": "code", "execution_count": 15, "metadata": {}, "outputs": [], "source": [ "wafer_thickness = 180 # um\n", "s_rear = 46 # cm/s\n", "cell_area = 243 # cm^2" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ " Other device parameters \n", "Other required device parameters: base diffusivity (in cm2/s), and optical generation profile, which allow us to estimate current collection in the device." ] }, { "cell_type": "code", "execution_count": 16, "metadata": {}, "outputs": [], "source": [ "generation_df = pd.read_excel(os.path.join(DATA_DIR, 'PVL_GenProfile.xlsx'), header = 0) # this is an optical generation profile generated by PVLighthouse's OPAL2 default model for 1-sun, normal incident AM1.5 sunlight on a 180-um thick SiNx-coated, pyramid-textured wafer.\n", "generation = generation_df['Generation (cm-3s-1)']\n", "depth = generation_df['Depth (um)']\n", "\n", "d_base = 27 # cm^2/s electron diffusivity. See https://www2.pvlighthouse.com.au/calculators/mobility%20calculator/mobility%20calculator.aspx for details" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Degradation parameters\n", "To model the device's degradation, we need to define several more important quantities about the degradation the device will experience. These include undegraded and degraded lifetime (in $\\mu s$)." ] }, { "cell_type": "code", "execution_count": 17, "metadata": {}, "outputs": [], "source": [ "tau_0 = 115 # us, carrier lifetime in non-degraded states, e.g. LETID/LID states A or C\n", "tau_deg = 55 # us, carrier lifetime in fully-degraded state, e.g. LETID/LID state B" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Remaining degradation parameters: \n", "\n", "The rest of the quantities to define are: the initial percentage of defects in each state (A, B, and C), and the dictionary of mechanism parameters.\n", "\n", "In this example, we'll assume the device starts in the fully-undegraded state (100% state A), and we'll use the kinetic parameters for LETID degradation from Repins." ] }, { "cell_type": "code", "execution_count": 18, "metadata": {}, "outputs": [], "source": [ "# starting defect state percentages\n", "nA_0 = 100\n", "nB_0 = 0\n", "nC_0 = 0\n", "\n", "mechanism_params = utilities.get_kinetics('repins')\n", "\n", "timesteps[['NA', 'NB', 'NC', 'tau']] = np.nan # create columns for defect state percentages and lifetime, fill with NaNs for now, to fill iteratively below\n", "\n", "timesteps.loc[0, ['NA', 'NB', 'NC']] = nA_0, nB_0, nC_0 # assign first timestep defect state percentages\n", "timesteps.loc[0, 'tau'] = letid.tau_now(tau_0, tau_deg, nB_0) # calculate tau for the first timestep" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Run through timesteps\n", "Since each timestep depends on the preceding timestep, we need to calculate in a loop. This will take a few minutes depending on the length of the timeseries." ] }, { "cell_type": "code", "execution_count": 19, "metadata": {}, "outputs": [], "source": [ "for index, timestep in timesteps.iterrows():\n", "\n", " # first row tau has already been assigned\n", " if index == 0:\n", " #calc device parameters for first row\n", " tau = tau_0\n", " jsc = collection.calculate_jsc_from_tau_cp(tau, wafer_thickness, d_base, s_rear, generation, depth)\n", " voc = letid.calc_voc_from_tau(tau, wafer_thickness, s_rear, jsc, temperature = 25)\n", " timesteps.at[index, 'Jsc'] = jsc\n", " timesteps.at[index, 'Voc'] = voc\n", "\n", " # loop through rows, new tau calculated based on previous NB. Reaction proceeds based on new tau.\n", " else:\n", " n_A = timesteps.at[index-1, 'NA']\n", " n_B = timesteps.at[index-1, 'NB']\n", " n_C = timesteps.at[index-1, 'NC']\n", "\n", " tau = letid.tau_now(tau_0, tau_deg, n_B)\n", " jsc = collection.calculate_jsc_from_tau_cp(tau, wafer_thickness, d_base, s_rear, generation, depth)\n", "\n", " temperature = timesteps.at[index, 'Temperature']\n", " injection = timesteps.at[index, 'Injection']\n", "\n", " # calculate defect reaction kinetics: reaction constant and carrier concentration factor.\n", " k_AB = letid.k_ij(mechanism_params['v_ab'], mechanism_params['ea_ab'], temperature)\n", " k_BA = letid.k_ij(mechanism_params['v_ba'], mechanism_params['ea_ba'], temperature)\n", " k_BC = letid.k_ij(mechanism_params['v_bc'], mechanism_params['ea_bc'], temperature)\n", " k_CB = letid.k_ij(mechanism_params['v_cb'], mechanism_params['ea_cb'], temperature)\n", "\n", " x_ab = letid.carrier_factor(tau, 'ab', temperature, injection, jsc, wafer_thickness, s_rear, mechanism_params)\n", " x_ba = letid.carrier_factor(tau, 'ba', temperature, injection, jsc, wafer_thickness, s_rear, mechanism_params)\n", " x_bc = letid.carrier_factor(tau, 'bc', temperature, injection, jsc, wafer_thickness, s_rear, mechanism_params)\n", "\n", " # calculate the instantaneous change in NA, NB, and NC\n", " dN_Adt = (k_BA * n_B * x_ba) - (k_AB * n_A * x_ab)\n", " dN_Bdt = (k_AB * n_A * x_ab) + (k_CB * n_C) - ((k_BA * x_ba + k_BC * x_bc) * n_B)\n", " dN_Cdt = (k_BC * n_B * x_bc) - (k_CB * n_C)\n", "\n", " t_step = (timesteps.at[index, 'Datetime'] - timesteps.at[index-1,'Datetime']).total_seconds()\n", "\n", " # assign new defect state percentages\n", " timesteps.at[index, 'NA'] = n_A + dN_Adt*t_step\n", " timesteps.at[index, 'NB'] = n_B + dN_Bdt*t_step\n", " timesteps.at[index, 'NC'] = n_C + dN_Cdt*t_step\n", "\n", " #calculate device parameters\n", " timesteps.at[index, 'tau'] = tau\n", " timesteps.at[index, 'Jsc'] = jsc\n", " timesteps.at[index, 'Voc'] = letid.calc_voc_from_tau(tau, wafer_thickness, s_rear, jsc, temperature = 25)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Finish calculating degraded device parameters.\n", "Now that we have calculated defect states, we can calculate all the quantities that depend on defect states." ] }, { "cell_type": "code", "execution_count": 20, "metadata": {}, "outputs": [], "source": [ "timesteps['tau'] = letid.tau_now(tau_0, tau_deg, timesteps['NB'])\n", "\n", "# calculate device Jsc for every timestep. Unfortunately this requires an integration so I think we have to run through a loop. Device Jsc allows calculation of device Voc.\n", "for index, timestep in timesteps.iterrows():\n", " jsc_now = collection.calculate_jsc_from_tau_cp(timesteps.at[index, 'tau'], wafer_thickness, d_base, s_rear, generation, depth)\n", " timesteps.at[index, 'Jsc'] = jsc_now\n", " timesteps.at[index, 'Voc'] = letid.calc_voc_from_tau(timesteps.at[index, 'tau'], wafer_thickness, s_rear, jsc_now, temperature = 25)\n" ] }, { "cell_type": "code", "execution_count": 21, "metadata": {}, "outputs": [ { "data": { "text/html": [ "
\n", "\n", "\n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", "
DatetimeTemperatureInjectionNANBNCtauJscVocIscFFPmpPmp_norm
01999-01-01 08:30:00-07:003.3487740.006239100.0000000.0000000.000000e+00115.00000041.5909970.66632710.1066120.8409875.6634671.000000
11999-01-01 09:30:00-07:0012.3784390.03160799.9975470.0024530.000000e+00114.99692241.5909870.66632710.1066100.8409875.6634600.999999
21999-01-01 10:30:00-07:0016.5950810.04422099.9919530.0080476.072523e-09114.98990541.5909660.66632610.1066050.8409865.6634460.999996
31999-01-01 11:30:00-07:0017.4572400.04320199.9859290.0140712.755257e-08114.98235041.5909420.66632410.1065990.8409865.6634300.999994
41999-01-01 12:30:00-07:006.3927830.00712499.9856640.0143362.864303e-08114.98201841.5909410.66632410.1065990.8409865.6634290.999993
..........................................
42961999-12-31 12:30:00-07:009.5155770.00115142.02957553.7666324.203792e+0072.48454641.3858980.65692910.0567730.8393175.5450200.979086
42971999-12-31 13:30:00-07:0028.1197960.04449942.03023153.7656544.204115e+0072.48503441.3859010.65692910.0567740.8393185.5450220.979086
42981999-12-31 14:30:00-07:0023.3146720.03513042.03631153.7594334.204257e+0072.48813541.3859240.65693010.0567800.8393185.5450340.979088
42991999-12-31 15:30:00-07:0017.8905280.02789142.04887353.7468134.204313e+0072.49442541.3859720.65693210.0567910.8393185.5450580.979093
43001999-12-31 16:30:00-07:004.5523650.00118842.04944553.7462424.204314e+0072.49471041.3859740.65693210.0567920.8393185.5450590.979093
\n", "

4301 rows × 13 columns

\n", "
" ], "text/plain": [ " Datetime Temperature Injection NA NB \\\n", "0 1999-01-01 08:30:00-07:00 3.348774 0.006239 100.000000 0.000000 \n", "1 1999-01-01 09:30:00-07:00 12.378439 0.031607 99.997547 0.002453 \n", "2 1999-01-01 10:30:00-07:00 16.595081 0.044220 99.991953 0.008047 \n", "3 1999-01-01 11:30:00-07:00 17.457240 0.043201 99.985929 0.014071 \n", "4 1999-01-01 12:30:00-07:00 6.392783 0.007124 99.985664 0.014336 \n", "... ... ... ... ... ... \n", "4296 1999-12-31 12:30:00-07:00 9.515577 0.001151 42.029575 53.766632 \n", "4297 1999-12-31 13:30:00-07:00 28.119796 0.044499 42.030231 53.765654 \n", "4298 1999-12-31 14:30:00-07:00 23.314672 0.035130 42.036311 53.759433 \n", "4299 1999-12-31 15:30:00-07:00 17.890528 0.027891 42.048873 53.746813 \n", "4300 1999-12-31 16:30:00-07:00 4.552365 0.001188 42.049445 53.746242 \n", "\n", " NC tau Jsc Voc Isc FF \\\n", "0 0.000000e+00 115.000000 41.590997 0.666327 10.106612 0.840987 \n", "1 0.000000e+00 114.996922 41.590987 0.666327 10.106610 0.840987 \n", "2 6.072523e-09 114.989905 41.590966 0.666326 10.106605 0.840986 \n", "3 2.755257e-08 114.982350 41.590942 0.666324 10.106599 0.840986 \n", "4 2.864303e-08 114.982018 41.590941 0.666324 10.106599 0.840986 \n", "... ... ... ... ... ... ... \n", "4296 4.203792e+00 72.484546 41.385898 0.656929 10.056773 0.839317 \n", "4297 4.204115e+00 72.485034 41.385901 0.656929 10.056774 0.839318 \n", "4298 4.204257e+00 72.488135 41.385924 0.656930 10.056780 0.839318 \n", "4299 4.204313e+00 72.494425 41.385972 0.656932 10.056791 0.839318 \n", "4300 4.204314e+00 72.494710 41.385974 0.656932 10.056792 0.839318 \n", "\n", " Pmp Pmp_norm \n", "0 5.663467 1.000000 \n", "1 5.663460 0.999999 \n", "2 5.663446 0.999996 \n", "3 5.663430 0.999994 \n", "4 5.663429 0.999993 \n", "... ... ... \n", "4296 5.545020 0.979086 \n", "4297 5.545022 0.979086 \n", "4298 5.545034 0.979088 \n", "4299 5.545058 0.979093 \n", "4300 5.545059 0.979093 \n", "\n", "[4301 rows x 13 columns]" ] }, "execution_count": 21, "metadata": {}, "output_type": "execute_result" } ], "source": [ "timesteps = letid.calc_device_params(timesteps, cell_area = 243) # this function quickly calculates the rest of the device parameters: Isc, FF, max power, and normalized max power\n", "\n", "timesteps\n" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Note of course that all these calculated device parameters are modeled STC device parameters, not the instantaneous, weather-dependent values. This isn't a robust performance model of a degraded module." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Plot the results" ] }, { "cell_type": "code", "execution_count": 22, "metadata": {}, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAoMAAAHcCAYAAACtXqAxAAAAOXRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjguMiwgaHR0cHM6Ly9tYXRwbG90bGliLm9yZy8g+/7EAAAACXBIWXMAAA9hAAAPYQGoP6dpAADoz0lEQVR4nOzdd3iN5xvA8e9JZE8kkSGSiBmzZu0VYis6UGpX1Vb8bLVHS6s2VdSqqlXU3ipWEiL2iBGEGNmy398fp46eJsg+GffnunI57/M+53nvc4Rzn+d9hkpRFAUhhBBCCJEv6ek6ACGEEEIIoTuSDAohhBBC5GOSDAohhBBC5GOSDAohhBBC5GOSDAohhBBC5GOSDAohhBBC5GOSDAohhBBC5GOSDAohhBBC5GOSDAohhBBC5GOSDAohdObbb79FpVLpOgwhhMjXJBkUIh+6fPkyXbt2xcnJCSMjIxwdHfn888+5fPlyhtqdMWMG27dvz5wgdeR1gvrs2bO31jl69CgqleqtP7/99pumnff9NGzYEIAePXpgbm6udZ2GDRtq6unp6WFpaUnp0qXp1q0bBw4cyMq3QQiRjxTQdQBCiOy1detWOnfuTKFChejduzdubm7cvXuXlStX8scff/Dbb7/Rvn37dLU9Y8YMPv74Yz766KPMDTqHGjx4MNWrV09WXqtWLTw8PChRooSmLDIykv79+9O+fXs6dOigKS9SpMg7r1G0aFFmzpwJQFRUFLdu3WLr1q2sW7eOTz/9lHXr1mFgYJBJr0gIkR9JMihEPnL79m26detG8eLFOX78OLa2tppzQ4YMoV69enTr1g1/f3+KFy+uw0izVnR0NKamphlup169enz88cdvPV+xYkXN42fPntG/f38qVqxI165dU30NKyurZPVnzZrF4MGDWbx4Ma6ursyePTvtwQshxD/kNrEQ+ch3331HdHQ0y5cv10oEAWxsbFi2bBlRUVHMmTNHU96jRw9cXV2TtfXf8X4qlYqoqCjWrFmjubXZo0cPzfmTJ09SvXp1jI2NcXd3Z9myZSnGmJCQwNSpU3F3d8fIyAhXV1fGjh1LbGxssrqLFy+mXLlymlvdAwYMIDQ0VKtOw4YNKV++PD4+PtSvXx9TU1PGjh2bincr59LX1+enn37Cw8ODhQsXEhYWpuuQhBC5mPQMCpGP7Ny5E1dXV+rVq5fi+fr16+Pq6sru3bvT3PbatWvp06cPNWrU4MsvvwTA3d0dgEuXLtGsWTNsbW359ttvSUhIYNKkSSneIu3Tpw9r1qzh448/5ptvvuHMmTPMnDmTq1evsm3bNk29b7/9lsmTJ+Pp6Un//v25fv06S5Ys4dy5c/z9999at06fP39OixYt6NSpE127dn3vrdnUioiISHFsYeHChbN8Yoy+vj6dO3dmwoQJnDx5klatWmXp9YQQeZckg0LkE2FhYTx69Ih27dq9s17FihX5888/iYiIwMLCItXtd+3ala+++orixYsnu605ceJEFEXhxIkTFCtWDICOHTtSoUIFrXoXL15kzZo19OnThxUrVgDw9ddfY2dnx/fff8+RI0do1KgRISEhzJw5k2bNmrFnzx709NQ3OcqUKcPAgQNZt24dPXv21LQbHBzM0qVL6devX6pfT2r06tUrxfLHjx9jb2+fqddKSfny5QH17X8hhEgvuU0sRD4REREB8N4E7/X58PDwTLluYmIi+/bt46OPPtIkggBly5bFy8tLq+5ff/0FwPDhw7XKv/nmGwBNj+XBgweJi4tj6NChmkQQoG/fvlhaWibr2TQyMtJKDjPLxIkTOXDgQLKfQoUKZfq1UvJ69vHrv1shhEgP6RkUIp94neS9L3FIbdKYWiEhIbx69YqSJUsmO1e6dGlNAghw79499PT0tGbhAtjb22Ntbc29e/c09V4//98MDQ0pXry45vxrTk5OGBoaZsrr+bcKFSrg6emZ6e2mVmRkJJB5f1dCiPxJegaFyCesrKxwcHDA39//nfX8/f1xcnLC0tIS4K1j3xITEzM9xtcye7ydiYlJpraXUwQEBAAkS56FECItJBkUIh9p3bo1gYGBnDx5MsXzJ06c4O7du7Ru3VpTVrBgwWQzdIFkvW+QchJna2uLiYkJN2/eTHbu+vXrWscuLi4kJSUlq/vkyRNCQ0NxcXHR1Evp+XFxcQQGBmrO52WJiYls2LABU1NT6tatq+twhBC5mCSDQuQjI0eOxMTEhH79+vH8+XOtcy9evOCrr77C1NSUkSNHasrd3d0JCwvT6lF8/Pix1sze18zMzJIljvr6+nh5ebF9+3bu37+vKb969Sr79u3TqtuyZUsAfvzxR63yefPmAWhmzHp6emJoaMhPP/2EoiiaeitXriQsLCzPz6xNTExk8ODBXL16lcGDB2t6cYUQIj1kzKAQ+UjJkiVZs2YNn3/+ORUqVEi2A8mzZ8/YuHGjZkkYgE6dOvG///2P9u3bM3jwYKKjo1myZAmlSpXC19dXq/2qVaty8OBB5s2bh6OjI25ubtSsWZPJkyezd+9e6tWrx9dff01CQgILFiygXLlyWklmpUqV6N69O8uXLyc0NJQGDRpw9uxZ1qxZw0cffUSjRo0AdW/jmDFjmDx5Ms2bN6dt27Zcv36dxYsXU7169TQt6vw28+bNS7YwtZ6entYahSdOnCAmJibZcytWrKi14HRGhIWFsW7dOkC9WPbrHUhu375Np06dmDp1aqZcRwiRjylCiHzH399f6dy5s+Lg4KAYGBgo9vb2SufOnZVLly6lWH///v1K+fLlFUNDQ6V06dLKunXrlEmTJin//S/k2rVrSv369RUTExMFULp37645d+zYMaVq1aqKoaGhUrx4cWXp0qUpthEfH69MnjxZcXNzUwwMDBRnZ2dlzJgxSkxMTLK4Fi5cqJQpU0YxMDBQihQpovTv3195+fKlVp0GDRoo5cqVS/V78zqmlH709fUVRVGUI0eOvLUOoEyaNClZuyEhIW89pyiK0r17d8XMzCxZ7P9u19zcXClZsqTStWtXZf/+/al+TUII8S4qRfnXPRYhhBBCCJGvyJhBIYQQQoh8TJJBIYQQQoh8TJJBIYQQQoh8TJJBIYQQQoh8TJJBIYQQQoh8TJJBIYQQQoh8TJJBIUSqrV69GpVKxd27d99b19XVlR49emR5TEIIITJGkkEh8qjAwEAGDhxIqVKlMDU1xdTUFA8PDwYMGKC160d+kpiYyKpVq2jYsCGFChXCyMgIV1dXevbsyfnz55PVv3z5Ml27dsXJyQkjIyMcHR35/PPPuXz5cobicHV11dr/+d/Onz+PSqVi9erVWuUnT56kRYsWODk5YWxsTLFixWjTpg0bNmx463Vq1KiBSqViyZIlGYpXCJG3STIoRB60a9cuypcvz9q1a/H09OSHH35g/vz5tGjRgr/++ovKlStz7949XYeZrV69ekXr1q3p1asXiqIwduxYlixZwhdffIG3tzc1atQgKChIU3/r1q1UqVKFQ4cO0bNnTxYvXkzv3r05cuQIVapUSXFv5qyyefNm6tevz5MnTxgyZAgLFiyga9euvHz5khUrVqT4nJs3b3Lu3DlcXV1Zv359tsUqhMh9ZG9iIfKY13vWuri4cOjQIRwcHLTOz549m8WLF6Onl7++C44cOZK9e/fyww8/MHToUK1zkyZN4ocfftAc3759m27dulG8eHGOHz+Ora2t5tyQIUOoV68e3bp1w9/fn+LFi2d57N9++y0eHh6cPn0aQ0NDrXNPnz5N8Tnr1q3Dzs6OuXPn8vHHH3P37l1cXV2zPFYhRO6Tvz4NhMgH5syZQ1RUFKtWrUqWCAIUKFCAwYMH4+zsrFV++PBh6tWrh5mZGdbW1rRr146rV6++93qKojBt2jSKFi2KqakpjRo1eutt1NDQUIYOHYqzszNGRkaUKFGC2bNnk5SUpKlz9+5dVCoV33//PcuXL8fd3R0jIyOqV6/OuXPn0vhuqAUFBbFs2TKaNm2aLBEE0NfXZ8SIERQtWhSA7777jujoaJYvX66VCALY2NiwbNkyoqKimDNnTrriSavbt29TvXr1ZIkggJ2dXYrP2bBhAx9//DGtW7fGysrqnbeThRD5m/QMCpHH7Nq1ixIlSlCzZs1UP+fgwYO0aNGC4sWL8+233/Lq1SsWLFhAnTp18PX1fWeP0sSJE5k2bRotW7akZcuW+Pr60qxZM+Li4rTqRUdH06BBAx4+fEi/fv0oVqwYp06dYsyYMTx+/Jgff/xRq/6GDRuIiIigX79+qFQq5syZQ4cOHbhz5w4GBgZpeUvYs2cPCQkJdOvWLVX1d+7ciaurK/Xq1UvxfP369XF1dWX37t1piiO9XvfyBgUFaRLWdzlz5gy3bt1i1apVGBoa0qFDB9avX8/YsWOzIVohRK6jCCHyjLCwMAVQPvroo2TnXr58qYSEhGh+oqOjNecqV66s2NnZKc+fP9eUXbx4UdHT01O++OILTdmqVasUQAkMDFQURVGePn2qGBoaKq1atVKSkpI09caOHasASvfu3TVlU6dOVczMzJQbN25oxTV69GhFX19fuX//vqIoihIYGKgASuHChZUXL15o6u3YsUMBlJ07d6b5fRk2bJgCKH5+fu+tGxoaqgBKu3bt3lmvbdu2CqCEh4enOR4XFxelVatWKZ47d+6cAiirVq3SlK1cuVIBFENDQ6VRo0bKhAkTlBMnTiiJiYkptjFw4EDF2dlZ83eyf//+VL9+IUT+I7eJhchDwsPDATA3N092rmHDhtja2mp+Fi1aBMDjx4+5cOECPXr0oFChQpr6FStWpGnTpvz1119vvd7BgweJi4tj0KBBqFQqTXlKt2I3b95MvXr1KFiwIM+ePdP8eHp6kpiYyPHjx7Xqf/bZZxQsWFBz/LqX7s6dO6l4J7S9fl8sLCzeWzciIiJVdV+ff912VurVqxd79+6lYcOGnDx5kqlTp1KvXj1KlizJqVOntOomJCSwadMmPvvsM83fSePGjbGzs5OJJEKIFMltYiHykNcJSmRkZLJzy5YtIyIigidPntC1a1dN+etZxaVLl072nLJly7Jv3z6ioqIwMzNLdv71c0uWLKlVbmtrq5XIgXp2q7+/f7IxeK/9dyJEsWLFtI5ft/fy5csUn/8ulpaWwJtE711ev4fvq5vapDG9/p1cA3h5eeHl5UV0dDQ+Pj5s2rSJpUuX0rp1a65du6YZO7h//35CQkKoUaMGt27d0jy/UaNGbNy4kdmzZ+e7yUNCiHeTZFCIPMTKygoHBwcCAgKSnXs9hjA1C0ZnhaSkJJo2bcqoUaNSPF+qVCmtY319/RTrKYqS5muXKVMGgEuXLlG5cuV31n39Hr5vLUZ/f3+cnJw0iWZaGBsb8+rVqxTPRUdHa+qkxNTUlHr16lGvXj1sbGyYPHkye/bsoXv37gCa3r9PP/00xecfO3aMRo0apTlmIUTeJcmgEHlMq1at+Pnnnzl79iw1atR4b30XFxcArl+/nuzctWvXsLGxSbFX8N/PvXnzptYSKyEhIcl68Nzd3YmMjMTT0zPVryWztGjRAn19fdatW5eqSSStW7dmxYoVnDx5krp16yY7f+LECe7evUu/fv3SFY+LiwtXrlxJ8dzrv4fX7+27VKtWDVDf6geIiopix44dfPbZZ3z88cfJ6g8ePJj169dLMiiE0CL3CoTIY0aNGoWpqSm9evXiyZMnyc7/t2fNwcGBypUrs2bNGkJDQzXlAQEB7N+/n5YtW771Wp6enhgYGLBgwQKtdv87MxjUPVXe3t7s27cv2bnQ0FASEhJS8erSx9nZmb59+7J//34WLFiQ7HxSUhJz587VLDo9cuRITExM6NevH8+fP9eq++LFC7766itMTU0ZOXJkuuJp2bIlQUFBbN++Xas8NjaWn3/+GTs7O6pUqaIpP3ToUIrtvB7P+foW/7Zt24iKimLAgAF8/PHHyX5at27Nli1biI2NTVfcQoi8SXoGhchjSpYsyYYNG+jcuTOlS5fm888/p1KlSiiKQmBgIBs2bEBPT09riZLvvvuOFi1aUKtWLXr37q1ZWsbKyopvv/32rdeytbVlxIgRzJw5k9atW9OyZUv8/PzYs2cPNjY2WnVHjhzJn3/+SevWrenRowdVq1YlKiqKS5cu8ccff3D37t1kz3mfu3fv4ubmRvfu3ZNt3/Zfc+fO5fbt2wwePJitW7fSunVrChYsyP3799m8eTPXrl2jU6dOmvdwzZo1fP7551SoUIHevXvj5ubG3bt3WblyJc+ePWPjxo24u7trXUOlUtGgQQOOHj36zli+/PJLfvnlFz755BN69erFBx98wPPnz9m0aRMBAQH8+uuvWmsKtmvXDjc3N9q0aYO7uztRUVEcPHiQnTt3Ur16ddq0aQOobxEXLlyY2rVrp3jdtm3bsmLFCnbv3k2HDh3e8+4KIfIN3U5mFkJklVu3bin9+/dXSpQooRgbGysmJiZKmTJllK+++kq5cOFCsvoHDx5U6tSpo5iYmCiWlpZKmzZtlCtXrmjV+e/SMoqiKImJicrkyZMVBwcHxcTERGnYsKESEBCguLi4aC0toyiKEhERoYwZM0YpUaKEYmhoqNjY2Ci1a9dWvv/+eyUuLk5RlDdLy3z33XfJYgSUSZMmaY4vXbqkAMro0aNT9Z4kJCQoP//8s1KvXj3FyspKMTAwUFxcXJSePXumuOyKv7+/0rlzZ8XBwUExMDBQ7O3tlc6dOyuXLl1KVjciIkIBlE6dOqUqlpcvXyrDhg1T3NzcFAMDA8XS0lJp1KiRsmfPnmR1N27cqHTq1Elxd3dXTExMFGNjY8XDw0MZN26cZmmbJ0+eKAUKFFC6dev21mtGR0crpqamSvv27VMVoxAif1ApSjpGYwshRA6wePFiRo0axe3btylSpIhOY/nrr79o3bo1Fy9epEKFCjqNRQgh0kLGDAohcq0jR44wePBgnSeCr2Pp1KmTJIJCiFxHegaFEEIIIfIx6RkUQgghhMjHJBkUQgghhMjHJBkUQgghhMjHJBkUQgghhMjHZNFpICEhAT8/P4oUKSIbuAshhBC5RFJSEk+ePOGDDz6gQAFJadJL3jnAz88vVXu4CiGEECLnOXv2LNWrV9d1GLmWJIOgWaPs7NmzODg46DgaIYQQQqTG48ePqVGjRo5YazQ3k2QQNLeGHRwctPZrFUIIIUTOJ0O8MkbePSGEEEKIfEySQSGEEEKIfEySQSGEEEKIfEzGDKZSUlIScXFxug5D5wwMDNDX19d1GELke4mJicTHx+s6DCGylHzmZA9JBlMhLi6OwMBAkpKSdB1KjmBtbY29vT0qlUrXoQiR7yiKQnBwMKGhoboORYhsIZ85WU+SwfdQFIXHjx+jr6+Ps7Nzvp6xpCgK0dHRPH36FECW4RFCB14ngnZ2dpiamsoHpMiz5DMn+0gy+B4JCQlER0fj6OiIqamprsPRORMTEwCePn2KnZ2ddN8LkY0SExM1iWDhwoV1HY4QWU4+c7JH/u3mSqXExEQADA0NdRxJzvE6KZbxSkJkr9f/5uSLqchP5DMn60kymEpyK+YNeS+E0C35NyjyE/l9z3o6TQaPHz9OmzZtcHR0RKVSsX37dq3ziqIwceJEHBwcMDExwdPTk5s3b2rVefHiBZ9//jmWlpZYW1vTu3dvIiMjs/FVCCGEECI3eF/ekZKjR49SpUoVjIyMKFGiBKtXr05WZ9GiRbi6umJsbEzNmjU5e/as1vmYmBgGDBhA4cKFMTc3p2PHjjx58iSTXlXG6TQZjIqKolKlSixatCjF83PmzOGnn35i6dKlnDlzBjMzM7y8vIiJidHU+fzzz7l8+TIHDhxg165dHD9+nC+//DK7XoIQQgghcon35R3/FRgYSKtWrWjUqBEXLlxg6NCh9OnTh3379mnqbNq0ieHDhzNp0iR8fX2pVKkSXl5emokvAMOGDWPnzp1s3ryZY8eO8ejRIzp06JDpry/dlBwCULZt26Y5TkpKUuzt7ZXvvvtOUxYaGqoYGRkpGzduVBRFUa5cuaIAyrlz5zR19uzZo6hUKuXhw4epvvaDBw8UQHnw4EGyc69evVKuXLmivHr1Kh2vSvfq16+vAMqGDRu0yn/66SfFwcEhXW3m9vdEiNxK/u3pRoMGDZQhQ4a89Tg7rpmfvev3/l2f3+/z37wjJaNGjVLKlSunVfbZZ58pXl5emuMaNWooAwYM0BwnJiYqjo6OysyZMxVFUecuBgYGyubNmzV1rl69qgCKt7d3muPOCjl2zGBgYCDBwcF4enpqyqysrKhZsybe3t4AeHt7Y21tTbVq1TR1PD090dPT48yZM9ke839duHiJoEePiYuPJylJQVHUP9lFURT8/PxwcHBgy5YtWud8fHyoUqVKtsUihBA9evRApVIxa9YsrfLt27fnqnFhW7duZerUqTqNISQkhP79+1OsWDGMjIywt7fHy8uLv//+G1CPs3vXz7fffguolyoaNGgQxYsXx8jICGdnZ9q0acOhQ4feef0GDRpo2jIwMKB06dJs2LAhq192tvP29tbKQwC8vLw0eUhcXBw+Pj5adfT09PD09NTU8fHxIT4+XqtOmTJlKFasmKaOruXYpWWCg4MBKFKkiFZ5kSJFNOeCg4Oxs7PTOl+gQAEKFSqkqZOS2NhYYmNjNccRERGZFbaWHr16cdH3vDouAwOMjIyxtLKieIlSVKj0AS6ubnRo/xFuzg7o66nQy+T/DG/evElERASzZs1i5MiRREdHa2Zl+fr60r59+0y9nhBCvI+xsTGzZ8+mX79+FCxYMFPajIuLy9YVHwoVKpRt13qbjh07EhcXx5o1ayhevDhPnjzh0KFDPH/+HIDHjx9r6m7atImJEydy/fp1TZm5uTl3796lTp06WFtb891331GhQgXi4+PZt28fAwYM4Nq1ayle+3VHw8yZM+nRowfR0dFMmTKF7t27U6tWLdzc3LL2xacgIiKC8PBwzbGRkRFGRkYZbjc4ODjFPCQ8PJxXr17x8uVLEhMTU6zz+v0LDg7G0NAQa2vrZHXelatkpxzbM5iVZs6ciZWVlebHw8Mj06+hKIrW9nUJ8fFERUbw+GEQfx87zNKf5jJm+EAqlS3J94t/Zt+pCwS9iCYyJj7Teg99fHwwNjamT58+WFpasmfPHkA9kPXq1avSMyiEyHaenp7Y29szc+bMt9aJjY1l8ODB2NnZYWxsTN26dTl37pzmfMOGDRk4cCBDhw7FxsYGLy8vTfmgQYMYOnQoBQsWpEiRIqxYsYKoqCh69uyJhYUFJUqU0Pxf+NrevXupW7cu1tbWFC5cmNatW3P79u23xtewYUOGDh0KwN27d1PseWvYsCGg3sp05syZuLm5YWJiQqVKlfjjjz+02ouKiuKLL77A3NwcBwcH5s6d+873MDQ0lBMnTjB79mwaNWqEi4sLNWrUYMyYMbRt2xYAe3t7zY+VlRUqlUqrzNzcnK+//hqVSsXZs2fp2LEjpUqVoly5cgwfPpzTp0+/9fqvOxrq1q2Lvb09xYsXZ8yYMSQkJODv7//O2LOKh4eH1uf6u36/RHI5Nhm0t7cHSDbb5smTJ5pz9vb2WgM0Qb1I9IsXLzR1UjJmzBjCwsI0P1euXEl1XIqiEB2X8N6fV/GJnPfxIeRlGIEPHnHlxm18/QPYe+goM7+bR6u2HwEQE/OKsYP70bJuFbp81pGTl27jHxTGrScR3H0WRdCLaILDXhH2Ki7NSaKvry8VK1bE0NCQ9u3ba/4DunjxIgkJCcmSwU8//ZQ6deqk6RpCCJEW+vr6zJgxgwULFhAUFJRinVGjRrFlyxbWrFmDr68vJUqUwMvLixcvXmjqrFmzBkNDQ/7++2+WLl2qVW5jY8PZs2cZNGgQ/fv355NPPqF27dr4+vrSrFkzunXrRnR0tOY5UVFRDB8+nPPnz3Po0CH09PRo3759qrYgdXZ25vHjx5ofPz8/ChcuTP369QF158Ovv/7K0qVLuXz5MsOGDaNr164cO3ZM08bIkSM5duwYO3bsYP/+/Rw9ehRfX9+3XtPc3Bxzc3O2b9+udZcrLV68eMHevXsZMGAAZmZmyc7/txfr33x8fFCpVFSsWFFT9vrv8r89ZNnlypUrWp/rY8aMyZR27e3tU8xDLC0tMTExwcbGBn19/ffmKnFxccm2kPx3HV3LsbeJ3dzcsLe359ChQ1SuXBmA8PBwzpw5Q//+/QGoVasWoaGh+Pj4ULVqVQAOHz5MUlISNWvWfGvb/+0+/nfX8vu8ik/EY+K+91d8p1JQthSF4pyJueNDzP1LJMVEcOLgXlocPoDdJ5Mxca2c7Fnbv66NhbEBhgX0sDAugKWJwTtvLfv6+moSvg4dOtChQwdiY2Px9fXF1tYWZ2dnTV1/f38ePHjArVu3MvjahBDZSflnyy5dSO92eO3bt6dy5cpMmjSJlStXap2LiopiyZIlrF69mhYtWgCwYsUKDhw4wMqVKxk5ciQAJUuWZM6cOcnarlSpEuPHjwfUX/xnzZqFjY0Nffv2BWDixIksWbIEf39/PvzwQ0B9y/XffvnlF2xtbbly5Qrly5d/52vR19fXfKDHxMTw0UcfUatWLb799ltiY2OZMWMGBw8epFatWgAUL16ckydPsmzZMho0aEBkZCQrV65k3bp1NGnSBFAntEWLFn3rNQsUKMDq1avp27cvS5cupUqVKjRo0IBOnTppJWjvcuvWLRRFoUyZMqmq/2++vr64ublhaWkJwLVr1xg1ahSVK1emRo0aAMybN49vvvmGgIAAypUrB8CAAQNYvHgxL168YNWqVcydOxcbGxuio6NZtGgRzZo1S3Msr1lYWGjiyUy1atXir7/+0io7cOCA5u/T0NCQqlWrcujQIT766CNA3Rt86NAhBg4cCEDVqlUxMDDg0KFDmt+169evc//+fU07uqbTZDAyMlIr+QgMDOTChQsUKlSIYsWKMXToUKZNm0bJkiVxc3NjwoQJODo6at7wsmXL0rx5c80/iPj4eAYOHEinTp1wdHTU0atKPYtKXlhU8kJJiCP0xDrCz26FpESebhqPkXN59M0Lo2dojKGdG2YVmhKXmERUXAJRcfAyOg5DfT2KFjLF1EAflSr5wpy+vr507twZUN/WMDAwYN++fSlOHpk4cSIzZszgyy+/5PHjx7IHpBC5RHR0NObm5jq5dmRkZIq9Sqkxe/ZsGjduzIgRI7TKb9++TXx8vNZdCgMDA2rUqMHVq1c1Za87AP7r38mQvr4+hQsXpkKFCpqy1z1X/76rdPPmTSZOnMiZM2d49uyZpkfw/v37700G/61Xr15ERERw4MAB9PT0uHXrFtHR0TRt2lSrXlxcHB988IHm9cbFxWl1YBQqVIjSpUu/81odO3akVatWnDhxgtOnT7Nnzx7mzJnDzz//TI8ePd4ba0aGI/n6+nL37l3Mzc1JSEhApVLx2WefMXv2bPT01DccAwICqFChAtevX6dcuXLcv38fb29vihYtSsGCBQkICGDu3Ll06tSJ7du38+2332YoGUyt9+UdY8aM4eHDh/z6668AfPXVVyxcuJBRo0bRq1cvDh8+zO+//87u3bs1bQwfPpzu3btTrVo1atSowY8//qgZmgDqya+9e/dm+PDhFCpUCEtLSwYNGkStWrU0X0h0TafJ4Pnz52nUqJHmePjw4QB0796d1atXM2rUKKKiovjyyy8JDQ2lbt267N27F2NjY81z1q9fz8CBA2nSpAl6enp07NiRn376KctiNjHQ58oUryxouQ1+vj40rFeXhIQEYh8EaJ2NPPYLq2P7Ua9hIwo5FMPBxZ24xCTuhKgX2NZTqTAsoKdOClERdC+Q0NBQipYoS3BYDNamBrRt25YtW7Zw6dIlzTdugHPnzhEZGUmjRo0oU6YMV65ckWRQCJGl6tevj5eXF2PGjElV8vJfb0tCDQwMtI5fz3b99zGgdQu4TZs2uLi4sGLFChwdHUlKSqJ8+fJa477fZ9q0aezbt4+zZ89iYWEBoNkAYffu3Tg5OWnVz4zJDcbGxjRt2pSmTZsyYcIE+vTpw6RJk1L1fpYsWRKVSvXWSSLv4uvry8iRI+nTpw+mpqY4ODgk64wICAjg008/1UxamTp1Ku3bt9fMdg4ICGDw4MEAODk5abZ+zWrvyzseP37M/fv3Nefd3NzYvXs3w4YNY/78+RQtWpSff/5ZM04V4LPPPiMkJISJEycSHBxM5cqV2bt3r9Yt8x9++EGTo8TGxuLl5cXixYuz4RWnjk6TwYYNG77z24lKpWLKlClMmTLlrXUKFSqUrdPZVSoVpoZZ87bV+bAmoaGhXLx4kcuXLxMZGUlQUBDLly8nMjKSJYsWsGTRAgBatWnH/6Z9j1UhGwCSFIWY+Df/mM77+GBgaIhtsZI8jYjhaUQMTVq0YVC/3kRHRzNu3DhN3YkTJ2qWSXidDL6+XSGEyNlMTU11tutSRvdInjVrFpUrV9bqBXN3d9eMBXRxcQHUe9KeO3dOM2kjMz1//pzr16+zYsUK6tWrB8DJkyfT1MaWLVuYMmUKe/bswd3dXVPu4eGBkZER9+/fp0GDBik+193dHQMDA86cOUOxYsUAePnyJTdu3Hjrc97Gw8MjVTtqgPqz08vLi0WLFjF48OBkyXVoaGiK4wbv3LlDaGgoTZs2pUSJEim2rSgKd+/epXXr1vzwww/cuXOHO3fuULNmTcqXL4+iKFy7do3SpUuTmJjIkiVLaNmyZZpea3q9L+9IaXeRhg0b4ufn9852Bw4cqLktnBJjY2MWLVqU6sWus1uOHTOYX5mZmVG7dm1q166tKZswYQKrVq3ir7/+4uDBgwDs3rkDkhL4888/UYD4RIX4xCQUBRQUHt66QlmPctgXsuB5pHqAcZmqdYhPSCQuLo6y5dW3Uv7++2+OHTumuZ0cGRnJJ598kr0vWgiRbiqVKt23anWtQoUKfP7551p3c8zMzOjfvz8jR47U3LqbM2cO0dHR9O7dO9NjKFiwIIULF2b58uU4ODhw//59Ro8enernBwQE8MUXX/C///2PcuXKaZYKMTQ0pFChQowYMYJhw4aRlJRE3bp1CQsL4++//8bS0pLu3btjbm5O7969GTlyJIULF8bOzo5x48Zpbrem5Pnz53zyySf06tWLihUrYmFhwfnz55kzZw7t2rVLdeyLFi2iTp061KhRgylTplCxYkUSEhI4cOAAS5Ys0bot/9rrySP/Xt/3vwIDA3F2dqZs2bLcvn2bqVOnMmHCBLZt20aVKlUIDAwkNjaWWrVqYWhoiKenZ5rec5H5JBnMBaytrRk2bBjDhg0jMjKSLl26sHPnTnbv3o2joyMzZ87E1taW5s2bU6CA+q/0h+/fDKy2NTfi1tNIMDLi1FV193cCcPNpBCNHj2XLn39Ru05tLI0NOHv2LKNGjdLFyxRC5ENTpkxh06ZNWmWzZs0iKSmJbt26ERERQbVq1di3b1+mrUv4b3p6evz2228MHjyY8uXLU7p0aX766SfN0jDvc/78eaKjo5k2bRrTpk3TlDdo0ICjR48ydepUbG1tmTlzJnfu3MHa2poqVaowduxYTd3vvvuOyMhI2rRpg4WFBd988w1hYWFvvaa5uTk1a9bkhx9+0IyxdHZ2pm/fvlrtvk/x4sXx9fVl+vTpfPPNNzx+/BhbW1uqVq3KkiVLUnzO69ndVlZWb203ICCA8uXLY2RkRFRUFA8fPqRhw4ZMnTqVL774goCAANq0aZNsiR2hOyolO7fEyKGCgoJwdnbmwYMHyWZwxcTEEBgYiJubm9ZYRV2KjY2lYcOGKa4D5erqioWFBVZWVlhbW2Nra4uNjQ3du3enZOmyRMYmEPRSPfPwzN/HWfHT9/y86U8AHKxM0IuLpEyZMoSEhLz1+jnxPREiP5B/eyI3mDFjBgYGBowcOZJNmzZRrlw5ypcvj4ODA3fu3OGHH34gPj6eSZMmpaq9d/3ev+vzW6Se9AzmQkZGRnh7e7NgwQIOHz7Mw4cP8fHxISkpibt376b4nDVr1nDnzh0KmppiVMCcxCSFIi2a0bRJE8Jj4gF4HPYK0Mf7ciAx8YkYG+hn34sSQgiRJwQEBNCtWzdAPbkC1LO3LSwsMDExISAggI8//liXIYr/kGQwFxs0aBCDBg0CIDExkfv37xMcHExUVBTPnj3jzp077N+/n2PHjvH06VMaNmzIuXPnMDNS/7Vbmqhn2L2KTyQwJIqEf2bXRcclcONJBO625pq6QgghRGqkNKnTzs6OGzduvPW80C35pM8j9PX1cXNzS7Yn5OjRo/n444/Ztm0b58+fx9bWlkqVKvHLL79oZq6ZGOjj4WhJXEISIZGxmgknt0MicSlshpWJQbLrCSGEECJvyLHb0YnMoaenx9atWzVLJjx79oxDhw5pFtT8N8MCejhZm+Bm82Zm4r3nUTz7JzlUFCXT9k0WQgghRM4gyWA+sXXrVv7880/NDLkJEybQpk0bZs2axf3797UmjFgYG+Bu+2ZHg0ehr7j6OJwrj8K59DCMBy/e3FIWQgghRO4myWA+YWNjQ5s2bTT7egLs2rWLMWPG4OLigqOjI+fOndOcMzMqgIeDJQb66l+R+MQkEv/pFYyOSyQ6NntWixdCCCFE1pJkMJ9p0aIFO3fuZOXKlZrV/QESEhI02wS9VkBfjzL2FhS3NdfcPn6dHCbK7WIhhBAiT5AJJPmMSqWidevWgHpT9aSkJMaMGcOcOXMYNmwYO3bsYPbs2dSoUUNT39yoAOb/zCouYpnEg5AYomIT+X7fdYrZWXH5YThPI2KIjE3A3soEA30VthZGlLW35GJQKJcfhdOkjB39Gri/NS4hROrJ2F2Rn8jve9aTZDCf09PTo1mzZsyZo96x5OjRo9SsWZM1a9bwxRdfJKv/732Z9wQ85mFE0H9qvEzxOmcDX3D3eTSdqjtTsahVsk3NhRDvZ2CgntkfHR2NiYmJjqMRIntER6s3Snj9+y8ynySDgiZNmnDv3j1mz57N4sWLAejevTtNmjTByclJq66xgT6uhc0Ie6qHa2FTSjkaU8HJCkdrE6LjErkVEsmtp5E8DnvFgxevKO9kScDDcAA2nr3PxrP32dCnJrVL2GT76xQit9PX18fa2pqnT58CYGpqKl+sRJ6lKArR0dE8ffoUa2tr9PVlI4SsIsmgAKBYsWIsWrSILl26ULduXQDKlCnDo0ePsLCw0KprZKCPrYURK3vUeOeWWLEJiRjq6/EoLIbfzt5nweFbAHT5+Qx7htSjrINl1r0gIfIoe3t7AE1CKEReZ21trfm9F1lDksE8rkGDBhw/fpwNGzbQuXNnTfmCBQuYOXMmjx490qpfp04dBg0axIIFC4iMjGTHjh107do1Xdc2KqD+FudkbcI3zUrjaG3CmK2XAGgx/wT1S9kyvGkpKjtbp+/FCZEPqVQqHBwcsLOzIz4+XtfhCJGlDAwMpEcwG0gymIcpioKfnx8ODg5s2bJFKxn08fGhSpUqKT5vxowZbNu2jaCgILp160apUqU0E0oy4uOqRXkRFcd3+64DcPxGCMdvhDCnY0U++sCJBy+jKVrQRJNECiHeTl9fXz4khRCZQpaWycNu3rxJREQE48ePZ8+ePZpBuAC+vr5UrVo1xeeZm5szc+ZMzfHChQszJR4DfT0GNCrB7/1q0dSjiKZ81BZ/So3fQ5O5x+j7qw+xCbKGoRBCCJFdJBnMw3x8fDA2NqZPnz5YWlqyZ88eAGJiYrh69epbewYBunTpwpgxYwBYu3Yttra2VKhQAZVKxaRJk4iLi0t3XDXcCrHii2qs71MTfT3twe/Hb4RQevxeRmy+yFdrffj9/AOi4xLSfS0hhBBCvJskg+kVF/X2n/iYNNR9lbq66eDr60vFihUxNDSkffv2/PHHHwBcvHiRhIQETTI4b948ihYtSuXKlXFzc2PAgAHo6ekxYsQILC3VkzyePXtGQEAAAHPmzGHTpk3piunf6pSw4ea0FmzpX5uDwxtQsaiV5twfPkHsvRzMqD/88Zi4j29+v8iZO88zfE0hhBBCaJMxg+k1w/Ht50o2g883vzn+rgTER6dc16Uu9Nz95vjHChCdQtLzbViaQ/T19dUkfB06dKBDhw7Exsbi6+uLra0tzs7OAAQEBDB//nw6duzIq1evsLGxYfbs2RQqVIhHjx4xduxYvL29cXFx0SSUjx8/TnM8KdHTU1HVpSAAOwbUwfvOc45dDyEuMYmNZ+8TE6/eA3mLbxBbfINoV9mRmR0qaK13KIQQQoj0k0/UPMzX11czaaRhw4YYGBiwb9++ZJNHAgICGDduHAAXLlzA3d0dc3NzAMzMzJg/f76m7tSpU5k4cSJRUenrrXwXlUpFbXcbarur1yCc1KYc607fY8HhmzwJjwVgx4VHBIfFMMSzJI9CYyhWyJQaboUyPRYhhBAiv5BkML3GPnr7OdV/ZviNvPWOuv+5Uz/0Uvpj+pc7d+4QGhqqSfoKFChA27Zt2bJlC5cuXaJFixaAesbx9evX6dixI9HR0YSGhnLy5Mm3tvt6zcHIyMhMifN9un7oQtcPXbj6OJwW808AcCbwBV1WnNHU6fCBEzM7VpBZyEIIIUQ6SDKYXoZmuq/7Dj4+PhgaGlK+fHlNWceOHenWrRvR0dGansDAwEDKli3L6dOnAZg2bRo//PADS5YsSbFdGxt1r92ePXs4e/YsNWvW5NSpU1y+fJng4GBcXFxo2rQpxsbG2Nvbc/bsWY4cOUK9evWoXbt2ul9PWQdL/CY0ZeimCwS9jCY6LpHHYeqxmVv9HlKvlA3tPyia7vaFEEKI/EqSwTzK19eX8uXLY2hoqClr2rQpiYmJxMXFaXoMAwICKF26tKZOhQoV8PX1fWu7rq6umscNGjTA1taWkJCQ98ZjY2PD06dPM7R1VkEzQ9b0erPeYXhMPBW/3Q/A3WdvGZMphBBCiHeS2cR51MyZM/Hx8dEqMzIyIjw8HEVRcHNzA7STwcTERDZu3EiTJk3e2m7t2rWZMGGC5jgkJIQiRYrw4YcfamYep+TZs2e0atWK588zb0awpbEBQ5qUBGBvQLCsTyiEEEKkg/QM5nMBAQGcPHmSzZs3o1KpaNGiBf37939rfT09PcaOHUvBggWJiIjgww8/pHHjxhQooP5Vevr0KZcvX+bPP/9k06ZN9O3bl3379nHmzBn27NmDjY0NAQEBlCtXLlPiL2Kp3hv5+pMISo/fy6+9alC/lG2mtC2EEELkBypFURRdB6FrQUFBODs78+DBA4oW1R53FhMTQ2BgIG5ubhgbG+sowpwlre/JgwcP6N+/P7t3q5fQmT9/PoMHD86UWEIiYvly7Xn87odqylZ8UU1rhxMhhBB507s+v0XqyW1ikeWcnZ3ZtWsXw4cPB+B///sfYWFpXzcxJbYWRmz7ug7r+9TUlPX99Tx1Zh1m7LZLPAp99Y5nCyGEEEKSQZFtSpUqBah7FqdPn56pbdcpYcOW/rWwMVdPmHkY+ooNZ+7T99fzxCcmZeq1hBBCiLxEkkGRbT7++GPN4++++w4vLy9evHiRae1XdSnEiVGNWd6tKvVKqpfAufwonOG/X8y0awghhBB5jSSDItsULlyYy5cva473799P2bJlefUq827lmhjq06ycPT93r0ZxW/Wajd63n2Va+0IIIUReI8mgyFYeHh789ttv6Ompf/WePn3Kl19+SXx8PFevXs20xNCogD47BtQB4FlkHCM2X0TmSgkhhBDJSTIost1nn31GYmIiDRs2BGDdunUYGhri4eHBgAEDMu06FsYGfFDMGoA/fIJYfPS2JIRCCCHEf0gyKHRm/vz5ycpWrVqFt7d3pl3jty8/1Nwu/m7fdWrNPMydkOzZV1kIIYTIDSQZFDpTsWJFwsPDmTJlCsuWLdOU165dm759+9KpUyf8/Pwy1JtnVECfpV2rYv/P4tTB4TGsP3M/w7ELIYQQeYUkg0KnLCwsmDBhAl9++SXjxo3TlP/8889s2rSJKlWqMGnSJBISEtJ9jVJFLDjxv0a0quAAwMqTgVwKypx1DoUQQojcTpJBkWNMmzaNnTt30qVLF5o3b64pnzp1KgYGBsyfPz/dvYQG+nr0a1Bcc9xm4UmeR8ZmOGYhhBAit5NkUOQorVu3Zv369ezZs4fz589rnRs6dCgrV65Md9sVi1ozu2MFzXGbBSdlQokQQoh8T5LBPK5BgwaoVCo2btyoVb5gwQIcHR11FFXqVK1ala1bt9KkSRNNWd++fdm+fXu62/ysejFGNS8NwKOwGEKj4zMaphBCCJGrSTKYhymKgp+fHw4ODmzZskXrnI+PD1WqVNFRZKnXvn17Dh48yE8//aQpW716NXFxcRw9ehQ/P780t/l1wxLYWhgBsPTY7UyLVQghhMiNJBnMw27evElERATjx49nz549REdHa875+vpStWpVHUaXNoMGDWLWrFkA7NixAyMjIxo1akT16tW5c+dOmttzs1EvN7Ps+B2uPg7P1FiFEEKI3ESSwXSKjo9O809C0psZsQlJCUTHRxOTEJOqdtPDx8cHY2Nj+vTpg6WlJXv27AEgJiaGq1ev5oqewX/r1q0bpUuX1ipLTEzE3d2dW7dupamtCa08NI9bzD/BJ0tPsejILZlUIoQQIt8poOsAcquaG2qm+TnfN/geL1cvAA7dP8SIYyOoVqQaq5qv0tRpvqU5L2NfJnvupe6X0nw9X19fKlasiKGhIe3bt+ePP/6gY8eOXLx4kYSEBK1kcPXq1SxYsIDExESSkpL45ptv6N69e5qvmZUcHR25du0aoaGhREREcPDgQXr16gXAoUOHKFGiRKrbqlDUirEtyzDjr2sAnLv7knN3X7Lz4iP2Dq2fJfELIYQQOZH0DOZhvr6+moSvQ4cO7N69m9jYWHx9fbG1tcXZ2RmAZcuW8csvv7Bv3z4uXLjAkSNHcvQsW2tra5ydnenZsyedOnUC4NmzZ2lu58v67hwf2YiRXqWxMTcE4FpwBDHxiZkarxBCCJGTSc9gOp3pcibNzzHUN9Q8blKsCWe6nEFPpZ2P7+24N8Oxvebr60vnzp0BaNiwIQYGBuzbt09r8sjLly8ZP348Fy9exMbGBoDChQvTo0ePTIsjK72eER0Wlr5FpIsVNmVAoxJ83dCdshP3EhOfxLXgCCo7W2dilEIIIUTOJclgOpkamGbo+QX0ClBAL/nbn9F2X7tz5w6hoaGapK9AgQK0bduWLVu2cOnSJVq0aAHAtm3baNSoUY5fZuZtrK2tgfQng6+pVCrsLY25+zyajxb9zblxnpoZx0IIIUReJreJ8ygfHx8MDQ0pX768pqxjx478+eefXL58WZMkBgQEULlyZR1FmXFWVlYALF++nKtXr2aora4fumget114kvAY9RqESUkKt55GcP7uC/ZcekyYrE0ohBAiD5GewTzK19eX8uXLY2j45tZ006ZNSUxMJC4uTpMMmpmZkZSUpKswM+z1uEcADw8P/v77b2rXrp2utvrUK86T8BhWnAjkcVgM9WYfoUV5ew5de0pIhPYs4+XdqtKsnH2GYhdCCCFyAukZzKNmzpyJj4+PVpmRkRHh4eEoioKbmxsALVq0YMOGDTx//hyA8PBw1q1bl+3xplebNm2YPXu25vj17OL0GtSkpGa8YNireH4794CQiFiMDfQwKvDmn8u03RnrhRRCCCFyCukZzOdq167N8OHDadSoEYqioFKpGDZsmK7DSrUCBQowatQonJ2d6dKlC9evX6djx478/vvv6Ovrp7k9S2MDfvvyQ344eIPwVwkYG+hR292GOiUKY2pYgJ0XHzFoox/3X0Tz58VHtK2UO8daCiGEEK+plJy8hkg2CQoKwtnZmQcPHlC0aFGtczExMQQGBuLm5oaxsbGOIsxZcup7UrFiRS5dUq/HWKJECWbOnEm7du0wMDDI1OvUmH6Qp//cNvab0JSCZobveYYQQois8K7Pb5F6cptY5Bl///03rq6uANy6dYtPPvmEH374IdOvs7DLm8W6P5h6gE+XetNv7Xk6Lffm0NUnmX49IYQQIitJMijyDAsLC06cOMFXX32lKdu9e3emX6eGWyF61HbVHJ+9+4J9l59w+s4Leq85z+FrkhAKIYTIPSQZFHlK0aJFWbJkCQcPHgTg+PHjjB07NtOvM7G1B3uH1uObpqVoU8kRDwdLzbleq89zPThCq35UbALet5/L7iZCCKFjixYtwtXVFWNjY2rWrMnZs2ffWjc+Pp4pU6bg7u6OsbExlSpVYu9e7c0hIiIiGDp0KC4uLpiYmFC7dm3OnTunVadHjx6oVCqtn+bNm2fJ60sPmUAi8qSaNd/sHb1t2zZmzJiRqe3r6akoY29JGfs3SeChq0/oveY8AJcfhVHa3oKY+ET+vPiIHw/c4FFYDAAOVsY0KWtHn7rFcbUxy9S4hBBCvN2mTZsYPnw4S5cupWbNmvz44494eXlx/fp17OzsktUfP34869atY8WKFZQpU4Z9+/bRvn17Tp06xQcffABAnz59CAgIYO3atTg6OrJu3To8PT25cuUKTk5OmraaN2/OqlWrNMdGRjlnYwOZQELqJpC4urpiYmKiowhzllevXnH37t0cN4Hkv65fv06ZMmWwtLTM8A4lqTVi80X+8AkCoMMHTlwICuVOSNRb69craUPLCg4kJClExyYQHZdITHwi8YkKHo6WtKrggIlh2mdFCyFEfpDWCSQ1a9akevXqLFy4EICkpCScnZ0ZNGgQo0ePTlbf0dGRcePGMWDAAE1Zx44dMTExYd26dbx69QoLCwt27NhBq1atNHWqVq1KixYtmDZtGqDuGQwNDWX79u0ZfMVZQ3oG3+P18iRxcXGSDP4jOjoaINNn6WY2e3v1otDh4eH4+vpqFtrOSjVcC2mSwa1+DwEwM9Snc41iDPEsSXhMArsuPmLmnmsAnLj5jBM3n721vb8uPeaXHtWzPG4hhMjr4uLi8PHxYcyYMZoyPT09PD098fb2TvE5sbGxyTo9TExMOHnyJAAJCQkkJia+s85rR48exc7OjoIFC9K4cWOmTZtG4cKFM+OlZViOTgYTExP59ttvWbduHcHBwTg6OtKjRw/Gjx+PSqUCQFEUJk2axIoVKwgNDaVOnTosWbKEkiVLZkoMBQoUwNTUlJCQEAwMDNDTy7/DLBVFITo6mqdPn2JtbZ2udfyyk6WlJaampkRHR/Phhx8SFRWV5Qnsp9WdcbUxY/QWf2zMjfD0sOOjD5yws1D/R2FhbEC/Bu60rezI7+eCOHL9KYb6elibGmBuVAATQ30MC+hx4MoTgl6+4vC1pzwOe4WpQQFeRsfxPCqW0Oh4XG3McLc1z9LXIoQQuUVERATh4eGaYyMjo2S3YZ89e0ZiYiJFihTRKi9SpAjXrl1LsV0vLy/mzZtH/fr1cXd359ChQ2zdupXERPX4bwsLC2rVqsXUqVMpW7YsRYoUYePGjXh7e1OiRAlNO82bN6dDhw64ublx+/Ztxo4dS4sWLfD29s4Rn6U5+jbxjBkzmDdvHmvWrKFcuXKcP3+enj17Mn36dAYPHgzA7NmzmTlzJmvWrMHNzY0JEyZw6dIlrly5kupbmO/rZo6LiyMwMDBXb9uWmaytrbG3t9ck5DnZr7/+Svfu3QHw9/enQoUKOo4odeISkig7cS+JSW//51ndtSCf13ShVUUHDPTz75cUIUT+9frz+78mTZrEt99+q1X26NEjnJycOHXqFLVq1dKUjxo1imPHjnHmzJlk7YSEhNC3b1927tyJSqXC3d0dT09PfvnlF169egXA7du36dWrF8ePH0dfX58qVapQqlQpfHx8uHo15d2q7ty5g7u7OwcPHqRJkyYZeAcyR45OBlu3bk2RIkVYuXKlpuzf9+oVRcHR0ZFvvvmGESNGABAWFkaRIkVYvXo1nTp1StV1UjPmICkpibi4uIy/qFzOwMAgR3yLSYuqVavi6+sLoBn/mRusPBnIr953ufdcfVve1FCfwuaGPHjxSquevp6KWsUL07+hO3VK2OgiVCGE0InXn9//nayRUs9gXFwcpqam/PHHH3z00Uea8u7duxMaGsqOHTveep2YmBieP3+Oo6Mjo0ePZteuXVy+fFmrTlRUFOHh4Tg4OPDZZ58RGRn5zuXNbG1tmTZtGv369Uvjq858Ofo2ce3atVm+fDk3btygVKlSXLx4kZMnTzJv3jxA/cEeHByMp6en5jlWVlbUrFkTb2/vtyaDsbGxxMbGao4jIiJSrPdvenp6OXqyhHi7Tz/9VJMMurm54e3tzYcffqjjqN6vd103etd1Iyw6HmNDPYwKqJPwsOh4fjt3n21+D7kWHEFiksLJW884eesZLoVN2TWoLhbGOXs8pxBCZCYLCwssLS3fWcfQ0JCqVaty6NAhTTKYlJTEoUOHGDhw4Dufa2xsjJOTE/Hx8WzZsoVPP/00WR0zMzPMzMx4+fIl+/btY86cOW9tLygoiOfPn+Pg4PD+F5cNUpUMDh8+PM0Njx8/nkKFCqX5ef82evRowsPDKVOmDPr6+iQmJjJ9+nQ+//xzAIKDgwFSvP//+lxKZs6cyeTJkzMUm8g9/ve//1GxYkVatmwJwB9//JErksHXrEwNkh33a+BOvwbuBDwMY/P5B2w6/4CY+CTuPY9m07kH9KlXXEfRCiFEzjV8+HC6d+9OtWrVqFGjBj/++CNRUVH07NkTgC+++AInJydmzpwJwJkzZ3j48CGVK1fm4cOHfPvttyQlJTFq1ChNm/v27UNRFEqXLs2tW7cYOXIkZcqU0bQZGRnJ5MmT6dixI/b29ty+fZtRo0ZRokQJvLy8sv9NSEGqksEff/yRWrVqYWiYuj1YT548ycCBAzOcDP7++++sX7+eDRs2UK5cOS5cuMDQoUNxdHTUjANLjzFjxmgluA8fPsTDwyNDsYqcrUWLFvzvf/9j9uzZPHz4UNfhZJryTlaUd7JiUpty9Fh9juM3Qpi2+yrFbc1oXKbI+xsQQoh85LPPPiMkJISJEycSHBxM5cqV2bt3r6ZT6f79+1oTRWNiYhg/fjx37tzB3Nycli1bsnbtWqytrTV1wsLCGDNmDEFBQRQqVIiOHTsyffp0zYRFfX19/P39WbNmDaGhoTg6OtKsWTOmTp2aY9YaTNWYQT09PYKDg1NckDElFhYWXLx4keLFM9Y74ezszOjRo7XW95k2bRrr1q3j2rVrmgGYfn5+VK5cWVOnQYMGVK5cmfnz56fqOrLRdf7w22+/0blzZwCuXLlC2bJldRxR5roUFEabhW+WMljdszoNS6fu36wQQuRG8vmdOVI1BXHVqlVYWVmlutFly5Ylu3WbHtHR0cmWctHX19fM6nVzc8Pe3p5Dhw5pzoeHh3PmzBmtmUJCgPr35bX3jQ/JjSoUtWLXoLqa4x6rztFt5Rk+XebNl7+e50WUTIASQgiRXKqSwe7du6epK7NLly6YmWV8m602bdowffp0du/ezd27d9m2bRvz5s2jffv2AKhUKoYOHcq0adP4888/uXTpEl988QWOjo5aM4WEAKhRowa9e/cG4PDhw9y8eVPHEWW+8k5WbB9QR3N84uYzzga+YP+VJ0zbdYWkdyxVI4QQIn/K0NIyAQEBHDt2jMTEROrUqUPVqlUzMzYiIiKYMGEC27Zt4+nTpzg6OtK5c2cmTpyoGb/4etHp5cuXExoaSt26dVm8eDGlSpVK9XWkmzn/iI6O1nxRsbe35/HjxzqOKGtcCgrj6uNwDAqo+PlEIJcfqRdjHd60FIObZM6C7EIIoWvy+Z050p0MLlq0iClTptCgQQPi4+M5fPgwo0aNYty4cZkdY5aTX6b8Ze7cuZp1KSMjIzOlFzsnexz2irqzj2gWsD44vAEl7GT3EiFE7ief35kj1cnggwcPtFb5Llu2LCdOnMDGRr3Irbe3N23btiUkJCRrIs1C8suUvyiKgomJCbGxsdy9excXFxddh5Tlbj2NxHPeMc1x3RI2rO5ZnQKyc4kQIheTz+/MkepPAk9PT+bPn8/r3LFw4cLs3buX2NhYIiIiOHjwILa2tlkWqBCZRaVSab7E5MYvL+lRws6c+Z0qU9hMPbzi5K1nDPv9It63n5OQKNssCiFEfpbqZPDcuXNcv36dmjVrcuHCBZYvX84PP/yAiYkJ1tbWbNq0iTVr1mRlrEJkmtfJYPXq1YmMjNRxNNmjXWUnzoxtQqki6lvEOy8+ovOK08w7cEPHkQkhhNClVG9HZ2lpyeLFizl16hQ9evSgcePGnDhxgsTERBITE7UWYBQip2vcuDEXL14E1F90GjVqpOOIskcBfT2WdK3KzycC2eX/iIiYBBYfvU1hcyN613V7fwNCCCHynDQPGKpduzbnz5+nYMGCfPDBBxw/flwSQZHrzJ07V7Po9ODBg7l69SoZmFifq7jbmjOzQwXOjXuzp/eGM/d0GJEQQghdSnUymJCQwOLFixk0aBCrV69m7Nix7Ny5k7lz5/LJJ5/w5MmTrIxTiEylUqmoW1e9QHNAQAAeHh44Oztz6tQpHUeWfYwN9Nk3tD4Aj8Ni8k0yLIQQQluqk8HevXuzcOFCzMzMWLVqFcOGDaNUqVIcPnyY5s2bU6tWLZYsWZKVsQqRqSZOnEibNm00xw8fPqRTp046jCj7uRQ2BSA6LpHN54N0HI0QQghdSHUyuGPHDrZs2cKsWbM4cOAAu3fv1pzr3bs3p0+f5sSJE1kSpBBZoWjRovz555/Ex8czb948QL2EUnx8vI4jyz7GBvoUsVTvLjRqiz8LDt3kUlCYjqMSQgiRnVKdDBYpUoT9+/cTFxfH4cOHKVy4sNZ5Ozs7NmzYkOkBCpHVChQowODBgzXHoaGhugtGB1Z2r655PPfADdosPMnm8w90GJEQQojslOpkcOHChUyfPh0TExO++uorfvzxxywMS4jspa+vj5WVFQAvXrzQcTTZq7yTFefGeTK4cQlN2e5Lj2UMoRBC5BOpTgabNm3KkydPCA4OJigoiNq1a2dlXEJku4IFCwLq2cX5LRGytTBieLPSLOzyAQBHr4cwZ991HUclhBAiO6RpaRmVSiW7jIg8y81Nvc7e/v37uXDhgm6D0ZF6JWwx/GeLuiVHb/PnxUc8j4zVcVRCCCGyUqqSwSpVqvDy5ctUN1q3bl0ePnyY7qCE0IW1a9dqHv/yyy/5rncQwMrUgJOjG6FSqY8Hb/Sj6rSD1J19mM7LT3Pq9jPdBiiEECLTpWoHkgsXLnDx4kUKFSqUqkYvXLhAbKz0JojcxcnJiZ49e7Jq1SoWLlyIp6cn7dq14+LFi6xevZobN25w4sQJ5s+fT8+ePXUdbpaxszBmYecq/BXwmKPXnhIVl0jQy1cEvXxFTEIi27620XWIQgghMpFKSUX3h56eHiqVKtU9JSqVips3b1K8ePEMB5gdgoKCcHZ25sGDBxQtWlTX4QgdOnfuHDVq1NAcN2vWjP379yert3DhQgYMGJCdoelEYpLCmcDn3AiO4NudVwCoX8qWRqVt6VHbFdXrLkQhhNAB+fzOHKlKBu/dS/tWVUWLFkVfXz9dQWU3+WUS/7Zu3Tq6deumVda4cWM++OAD5s6dC0Dx4sW5ffu2LsLTmeY/HudacITmePfgupRztNJhREKI/E4+vzNHqm4Tu7i4ZHUcQuQYH3/8Mf7+/rx48QJbW1s6depEpUqVAGjQoAFt27bl4cOHKIqSr3rGNn9Vi79vPeerdT4ArD9znxntK+g4KiGEEBmVqp7BvE6+WYjUio2NxdjYGFD3IH7++ec6jij7jd9+iXWn7wPg4WDJrI4VqFjUWrdBCSHyJfn8zhxpWlpGiPzOyMgIJycnALp27Up4eLiOI8p+PWq7UcFJfXv4yuNw+q/zJSY+UcdRCSGESC9JBoVIo3Xr1mke79y5U4eR6EYJO3N2DqrLvE/Vt84fhr6iydxjPJP1CIUQIleSZFCINGrYsCHVq6v38+3atSsmJiYULVqU8uXLc+LECR1Hl33af+BEhw/UvaQPQ19RZ9ZhImLidRyVEEKItEpzMvjgwQOCgoI0x2fPnmXo0KEsX748UwMTIifr3Lmz5nFMTAwPHz7k8uXLtGjRgpiYGB1Gln1UKhXzPqvMgs7qLexiE5K4FBSm46iEEEKkVZqTwS5dunDkyBEAgoODadq0KWfPnmXcuHFMmTIl0wMUIicaNmwY0dHRjB8/ngoV3syojYqK4tNPPyUpKUmH0WWvNpUcqVdSvRB1l5/PMGF7AC+j4nQclRBCiNRKczIYEBCgWZT3999/p3z58pw6dYr169ezevXqzI5PiBzLxMSEqVOn4u/vj6IoNGnSBFCPIzx16lSy+nm5x7D9P7eLAdaevse8Azd0GI0QQoi0SHMyGB8fj5GREQAHDx6kbdu2AJQpU4bHjx9nbnRC5CK//PKL5vG/h1IkJibSr18/zMzMqFOnDtOmTePWrVu6CDHLdKhSlIuTmlHJ2RpQJ4Sbzt3XbVBCZFTUM/BbD+dWQrh8vom8K83JYLly5Vi6dCknTpzgwIEDNG/eHIBHjx5RuHDhTA9QiNyiWLFifPzxxwA8e/YMgOvXr9OsWTOWL19OUlISp06dYsKECZQsWZL58+frMtxMZ2ViwJqe1TXHM/66RlJSvl/GVOjaq1A4Okv9c2a5OsFLiaJA+CN1nVuH1GXei2DH17B7OMwrA8sbwfP8tfOQyB9StQPJv82ePZv27dvz3Xff0b17d83ODH/++afWnq5C5Ec2Nuqxc4MGDUKlUjF27FjNWoSfffYZTk5OzJs3D4ChQ4dy8+ZNxo0bx4MHD3j69ClPnz7l2bNn6Onp0blzZ82ahrmFtakhR0c0pOH3Rwl7FU/lKfup5V6YmR0qUsjMUNfhifwm8DisaaNdVr7Dm8cv74GePgRshdNLIOKRurzDCvWfnpNApYIT6m0oeeSrTgj7HgabElkfvxDZJF07kCQmJhIeHk7BggU1ZXfv3sXU1BQ7O7tMDTA7yArmIrOsWrWKXr16aZVZWlqyYcMGWrVqBcDt27cpUSJ1HyRTp05l/PjxmR5nVhuz1Z+NZx9ojnvXdWNCaw8dRiTyjYQ4+HMQ3D8Fof8aqmDnAbal4ZPVb8qW1IUnl1Jup/cBcP6ngyP+FQRsAe/F8PSyuqzvEXCqkiUvQaRefvz8PnjwIJMmTSIpKYmWLVsyevRoDAwMMtRmutYZVBQFHx8fli1bRkSEeuN6Q0NDTE1NMxSMELldz549NROp9PX16dSpE6dPn9YkggDu7u5cu3aNhg0basqcnJyoWrUqLVq0oF27dpryFStWZFfomWpmh4pcnuxFnRLqoSMrTway+fyD9zxLiEwQFwkWRSDq+Zuyz9bB197aiWBsBDy7/ubYtqz6/NjH8G3Ym0QQwMAEPugKn28G08JQ0A3sK6rPHftOfQs6MQEe+qgT0KdX39xOjo+BQ1Pg2l/wIhCiX0DIdbh/5s0t6+e34cBEODZHfV6IdxgwYACjRo1iwYIF3L17l8mTJ2e4zTT3DN67d4/mzZtz//59YmNjuXHjBsWLF2fIkCHExsaydOnSDAeV3fLjNwuRtfz9/SlUqNB7f5+eP3+OiYlJsi9S169fp0yZMoC6J15P7833tqSkJFQqFSqVKvMDz2SPQl9Re9ZhQD3j+IfPKus2IJF3nVoI8dHQYJR6/N/902BaCKyLqZO5lCgKRASDhb36dnBqRD2HiMdgXx7uecOq5inX+2gJVO6iTu5+rAhxEcnrlOsAn6yCCxth+1fqMofK6tvQevqpiyefy4+f3x988AF+fn6A+vOhTp06nD59OkNtprlncMiQIVSrVo2XL19iYvLmH1j79u05dOhQhoIRIq+oWLFiqv5jKly4cIo96m5ubprHRkZGdOzYkXr16uHi4oK+vj729vY8eJDze9ocrU1Y1EV9K22b30N+ORmo44hEnhATDr5r4ehsdc/chk6wfxxEP4cXd9SJnUst9W3htyWCoK5n6ZD6RBDArLA6EQR1smleJOV6Fg7qP4MvvUkE9f8zbtapqvrPSp3A4587Ao8vwLmfUx+PyHdCQkLYvHkz/v7+JCQkEBeX8XVd09wzWLhwYU6dOkXp0qWxsLDg4sWLFC9enLt37+Lh4UF0dHSGg8pu+fGbhcj5OnTowLZt295Z55NPPmHx4sWaiSs50dOIGGpMV39RdC5kwolRjXUckcj1tvSBS5uTlxuYwjfXwdgye+N5chnioqGgi/o2saEZFHKHAv9K/uJjoICR+hZyVAgYW4FDJe1E9Mxy2DNS/bhKdyhaDSp8CgbG2ft6cpH8+Pk9b948Ll++TEBAANeuXSM2Npa2bdtSoUIFKlSowEcffZTmNtM8mzgpKYnExMRk5UFBQVhYWKQ5ACFEyrZu3Yqfnx+bNm3Czs6OokWL4uzszL59+zRjRDZv3szmzZvp1KkTY8aMoWLFijqOOjk7C2NOj2nChzMP8eDFK3ZefESbSo66DkvkVIoCkU/Ut24BQm7Aie/Vy7702KUusyur/rNoDSjsrr5ta1IQ6gzJ/kQQoEi5N4/N3zKJ8nVCV+QdE6mq94Ybe+D2YfBdA0mJUP7jzItT5Am9e/fGyspKcxwYGEhAQAABAQFs3rw5XclgmnsGP/vsM6ysrFi+fDkWFhb4+/tja2tLu3btKFasGKtWrUpzELqWH79ZiNwtJCSEdu3a4e3trVX+448/MmTIEB1F9W4NvjvCvefqOwfOhUz4qoE7n1R1xrBAuuaxidwuIhgOTIIHZ6BMK/D8FvQN1OsCznaFmv3gyg51ogdg7QJfnVQnew991DN8Xevq8AVkkYhg2DEAVHrQaq56zKN4q/z4+a2np4ezszPlypWjfPnyVKhQgfLly+Ph4aHZFCSt0pwMBgUF4eXlhaIo3Lx5k2rVqnHz5k1sbGw4fvy4LC0jRDbavXs3/fv31xo/ePHixRzZQ3j1cTgjNl/k8qNwTVkN10JMbOOBrYURhcwMMdCXxDBfeHYLfm4MMWFvymr2hxaz1Ang718kf06lLtBsmnrMXn6SEAvHv4d7f0PXrXLL+D/y4+f3gAEDOH36NB06dMDe3p5Lly4REBDAlStXsLS05Nq1a2luM13rDCYkJPDbb7/h7+9PZGQkVapU4fPPP9eaUJKb5MdfJpG3XL58mfLl1YPaN2zYQOfOnXUcUcoSEpM4dO0p03df5f4L7fHFdhZGNCtXhEJmRhgV0KOCkxU13AphbCCzKvOUxARYWBVe3lUfe7SD2kPAoaK6ZzD6hXpM4O1D6nF1npOhyhf5d3atosC8suoe0m7bwb2RriPKUfLr5/eDBw+YOnUqgYGBTJgwgfr16wPqu0a2trZpbi9dyWBek19/mUTe0qNHD9asWQPArVu3cHd313FEbxebkMikHZfxufeSkMhYQqPjU6xX2MyQw980xMo0YwuqCh1JTIBgf3UPoHsjdWLzS3N48M8yGHWGqm8PpzSbNyYMjCzTNtM3r9rWHy5uUL9fTTO+plxekt8/vwMDA5k2bRqPHj1i7ty5eHikb3H/NE8g+fPPP1MsV6lUGBsbU6JECa1lMYQQ2aNx48aaZPD3339nzJgxOo7o7YwK6DOr45tb2WGv4ln9913CXsUTn5jEjScRnAl8wfOoOH44eINJbTxyxbqK4l9uHYINn0JSAlTrrU4GkxKg9D/JYKnm705sjK3efi6/Kd5QnQwGbIVGY9WzkkW+dfXqVa5fv87169e5evUqt2/fJioqisuXL6c7GUxzz6Cenh4qlYr/Pu11mUqlom7dumzfvl1ru7qcLL9/sxB5x9ChQ5k/fz4Ae/bsoXnztyyImwv8e0u7pV2r0ry8vY4jEu8UcgOOzYKwh+revPv/mtxUa6B6vN/rhD4pCfRkfGiqxUbCgqoQGQwtv4cafXUdUY6RHz+/9fT0qFixIp9++imtW7embNmy2b8d3YEDB6hevToHDhwgLCyMsLAwDhw4QM2aNdm1axfHjx/n+fPnjBgxIkOBCSHSrkePHprHw4cPx9fXl5MnT3Lz5k3dBZVOXzV4c5t76bHbhL3lVrLQgbsn4WdPmFUMnv4zWD0+Cp5cUff6vU4EDUzhiz/Ba7r27V5JBNPGyBzq//OZum8cnF6ivuUu8qXvv/+eqlWrsmPHDpo2bUrt2rXp2bMnc+fOZd++felqM809g+XLl2f58uXUrl1bq/zvv//myy+/5PLlyxw8eJBevXpx//79t7SSs+THbxYi79qyZQsff5x8bbIjR45o7YecG2w8e58xWy8B0KO2K9+2LfeeZ4gslxAHMxwh6Z/kvO/hNztpBPmA9wL1ki9JSVDxE/X6fyLjEmJhYTX1otYydlBDPr+11xm8fPky69atS3MbaR4zePv2bSwtky/qaWlpyZ07dwAoWbIkz549S3MwQoiMa9u2LR9++CE+Pj5YWlry/PlzAE6dOpXrksGPqxZlrfc9rjwOZ/WpuzQsbUvD0rlv+ao8IykJFlV/kwg2mwZFKrw5X7QqfLJaJ6HleQWMoOdeuLwVan71ptxnNSTGQ/U+MtkmD7t37x7+/v4UKVKEGjVqaJ1zc3PDzc2NNm3apLv9NPfVV61alZEjRxISEqIpCwkJYdSoUVSvXh2Amzdv4uzsnO6ghBDpZ2BggLe3N7GxsTx79kyzW8m4cePYv3+/jqNLGwN9PVb2qKY57rHqHB8t+pvYhOS7IIksFvYQptu/WRKmel+oPUh7yzWRtayc1O+5/j/jw67vhV3D4K8R6t7avWPgwgZ4cE69RI/cSs4TNm7cSKlSpWjXrh21atWiWrVqWjlYZkhzMrhy5UoCAwMpWrQoJUqUoESJEhQtWpS7d+/y88/qzbUjIyMZP358pgYqhEib17NvmzZtqilbuXKlrsJJNwcrE/74qhZO1up1TC88CGX96dwxBCVPuf4XJMaqH9ceDK2+1208ApxrgG0Z9eP4aDi9GLb3h5WeMEdW9cgrJk+eTJcuXbh27ZrmC/3o0aMz9RrpWmcwKSmJ/fv3c+PGDQBKly5N06ZN0culg4JlzIHI61avXk3Pnj2xsLDg9u3b6VqUNCcY8psfOy48AmDzV7Wo7lpIxxHlI0mJ8PQKRDyBEk3klmROcuVP9Q4lAEHnISwIEuPgf4G6jSsb5IfPb0NDQ27cuIGrqysA165do2rVqkRFRWXaNdI8ZhDU05qbN2+eq5etECI/qVmzJgARERHY2dmxa9cuWrVqpeOo0m5wk5KaZPBSUJgkg9lJTx/sK6h/RM7i0Vb982/xr/71OEY95lAS+FwpISEBU1NTzXGZMmVISkoiODgYe/vMWXIrXclgVFQUx44d4/79+8TFxWmdGzx4cKYEJoTIPGXLlmX69OmMGzcOgDZt2vDkyZNc10PobmtOj9qurD51l5DIWF2Hkz/4rYNjc6DJRHBvDKaSgOcKBv9sD7t/PJxaAB1XQoXkqwyI3GHNmjXUqVOHihUrYm5uToECBYiOjn7/E1Mpzcmgn58fLVu2JDo6mqioKAoVKsSzZ88wNTXFzs5OkkEhcqixY8fy4Ycf0qRJExRFYfr06Xz//fcUKJCu74Q6Y2uh3n1hxfE79G/ojqWxbFWXZe55w44B6sdbesOwK7qNR6Rd5FP1n1v7QslmYJx8NRCRs9WrV49p06YRERGBnp4ebm5uxMTEsHLlSjw9PalWrRoWFhYZukaaB/kNGzaMNm3a8PLlS0xMTDh9+jT37t2jatWqfP+9DCgWIidr3LgxXl5eAMyfPx8DAwPWr1+v46jS5vVEkoQkhWXHbus4mjxIUeDWQVhUE1b9ayjQ12fUs1lF7tLwn4kGShJs6aOeZSxylWPHjhEWFsb169dZt24d7du3p0GDBixZsoQmTZpQsGBBypYtm6FrpDkZvHDhAt988w16enro6+sTGxuLs7Mzc+bMYezYsRkKRgiR9SZOnEjVqlU1xwsWLNBhNGnXvLw9RQuqE8KbTyJ1HE0edPE3WNcRQv7ZWcTADD5dC3ZldBuXSJ9CxdV/fwA396m3tTu9VL1mZD61aNEiXF1dMTY2pmbNmpw9e/atdePj45kyZQru7u4YGxtTqVIl9u7dq1UnIiKCoUOH4uLigomJCbVr1+bcuXNadRRFYeLEiTg4OGBiYoKnp2ead4YqWbIknTp1Ys6cORw8eJAXL15w+/ZtNm7cSPv27dPU1n+lORk0MDDQzBq2s7PT7DJiZWXFgwcPMhSMECLr1a5dm/Pnz2v+szpz5kyuWgrK2ECf8a3Um7E/k3GDma+IB5RuCYYW0GoejLqdfHKCyF3KtoEGo9XbA756AXv/B+dz3zJTmWHTpk0MHz6cSZMm4evrS6VKlfDy8uLp06cp1h8/fjzLli1jwYIFXLlyha+++or27dvj5+enqdOnTx8OHDjA2rVruXTpEs2aNcPT05OHDx9q6syZM4effvqJpUuXcubMGczMzPDy8iImJiZDr8fNzY1PPvmEGTNmZKgdlDRq2rSpsn79ekVRFKVPnz5KjRo1lHXr1ileXl5KjRo10tpcjvDgwQMFUB48eKDrUITINvHx8YqJiYkCKOXKldN1OGly/u5zxeV/u5Q6sw7pOpTcLz5GUf4crCi3j+g6EpHVol8qys6hinJ2haK8CNR1NJkirZ/fNWrUUAYMGKA5TkxMVBwdHZWZM2emWN/BwUFZuHChVlmHDh2Uzz//XFEURYmOjlb09fWVXbt2adWpUqWKMm7cOEVRFCUpKUmxt7dXvvvuO8350NBQxcjISNm4cWOq4s5qae4ZnDFjBg4ODgBMnz6dggUL0r9/f0JCQli2bFnGMlMhRLYpUKAAp0+fBsj01eyzmo25ehJJ0MtX/HQobbdaxH/cOabe0uzXdnB4mq6jEVnJxBpa/6Deuq6gq66jyXZxcXH4+Pjg6empKdPT08PT0xNvb+8UnxMbG4uxsbFWmYmJCSdPngTUy74kJia+s05gYCDBwcFa17WysqJmzZpvvW52S/M0wmrV3mwNZWdnl+zeuRAi9yhcuDAAL168QFEUza4lOZ2TtQlFC5oQ9PIV8w7coIilEZ9VL6brsHKf6Bew4RP1Y2sXqDVAt/EIkU4RERGEh4drjo2MjDAyMtKq8+zZMxITEylSpIhWeZEiRbh27VqK7Xp5eTFv3jzq16+Pu7s7hw4dYuvWrSQmqrfEtLCwoFatWkydOpWyZctSpEgRNm7ciLe3NyVKlAAgODhYc53/Xvf1OV1Lc89g48aNCQ0NTVYeHh5O48aNMyMmIUQ2KVRIvWZcQkICL17knlmGBfT1ODqiIZbG6u+zY7ZeIuxVvI6jykXioiHIB37r8qas9iAwKai7mET2CQuCsytgdWvw/13X0WQKDw8PrKysND8zZ87MlHbnz59PyZIlKVOmDIaGhgwcOJCePXtq7bi2du1aFEXByckJIyMjfvrpJzp37pyrdmVLc6RHjx5NttA0QExMDCdOnMiUoIQQ2cPExAQTE/XM3A8++AAlF21sX0Bfj61f1wEgSYEHLzJvAdY8b9dQ+Lkx3P/nFlWlLlCtl05DEtko/BHsGwd3T6jXHzyVu1YUSMmVK1cICwvT/IwZMyZZHRsbG/T19Xny5IlW+ZMnT966k4etrS3bt28nKiqKe/fuce3aNczNzSlevLimjru7O8eOHSMyMpIHDx5w9uxZ4uPjNXVet52W6/7b4cOH8fDw0Or5fC0sLIxy5cplOP9KdTLo7++Pv78/oH7TXx/7+/vj5+fHypUrcXLK/DWoHj58SNeuXSlcuDAmJiZUqFCB8+fPa84rmTBdW4j8rFcvdRLw4MEDPD09OX/+fK5JCkvYmVPeSb2IbnBYxmbl5QuKov5pvww8J4NTNXUi2Gquers5kT8414Bhl6FGP/Xx/vFw56hOQ8ooCwsLLC0tNT//vUUM6j1+q1atyqFDhzRlSUlJHDp0iFq1ar2zfWNjY5ycnEhISGDLli20a9cuWR0zMzMcHBx4+fIl+/bt09Rxc3PD3t5e67rh4eGcOXPmvdcF+PHHH+nbty+WlskXDLeysqJfv37Mmzfvve28U2pnmqhUKkVPT0/R09NTVCpVsh9TU1Nl5cqVmTq75cWLF4qLi4vSo0cP5cyZM8qdO3eUffv2Kbdu3dLUmTVrlmJlZaVs375duXjxotK2bVvFzc1NefXqVaqvI7OJRX7Xpk0bBdD8/P7777oOKdV6rz6nuPxvlzJgvY+uQ8m5nl5TlA2dFGW6k6LEhOs6GpGT/NFHUSZZKsq6T3QdSbqk9fP7t99+U4yMjJTVq1crV65cUb788kvF2tpaCQ4OVhRFUbp166aMHj1aU//06dPKli1blNu3byvHjx9XGjdurLi5uSkvX77U1Nm7d6+yZ88e5c6dO8r+/fuVSpUqKTVr1lTi4uI0dWbNmqVYW1srO3bsUPz9/ZV27dqlOlcpVqyYcuXKlbeev3r1quLs7Jyq1/82qZ5AEhgYiKIoFC9enLNnz2rtaWpoaIidnR36+pn7zXL27Nk4OzuzatUqTZmbm5vmsaIo/Pjjj4wfP16Tgf/6668UKVKE7du306lTp0yNR4i8as2aNXz77bf89NNPABw5coRPPvlEx1GljpO1ehbfLv/HdKr+jLolbXQcUQ60b6x6VxFQ3yK0La3beETO0XA0XPpdvSD19gHgNV096ziP+uyzzwgJCWHixIkEBwdTuXJl9u7dq5nccf/+fa2xfjExMYwfP547d+5gbm5Oy5YtWbt2LdbW1po6r29LBwUFUahQITp27Mj06dMxMHizVeaoUaOIioriyy+/JDQ0lLp167J3795ks5BT8uTJE622/qtAgQIZXhFCpSg5936Qh4cHXl5eBAUFcezYMZycnPj666/p27cvAHfu3MHd3R0/Pz8qV66seV6DBg2oXLky8+fPT7Hd2NhYYmPfLFb78OFDPDw8ePDgAUWLFs3S1yRETrZkyRK+/vprGjZsyJEjR3QdTqrcehqJ57xjADQvZ8+SrlVyzazobPHXKDj7z7JfdYZA44mgn7v2oxZZbP2n6mSw1kBoOhVy0cSHoKAgnJ2d8/Tnt7u7O3PnzuWjjz5K8fzWrVsZMWIEd+7cSfc10vU/ws2bNzly5AhPnz4l6T9b2kycODHdwfzXnTt3WLJkCcOHD2fs2LGcO3eOwYMHY2hoSPfu3dM9XXvmzJlMnjw50+IUIq9wdHQE1BPFpk+fzrhx43Qc0fuVsDNnYmsPpuy6wt7LwWw4e5/Pa7roOqycISkJ/P7Zisy8CDQaL4mgSK7dQri2Sz1+NBclgvlFy5YtmTBhAs2bN0/Wk/jq1SsmTZpE69atM3SNNPcMrlixgv79+2NjY4O9vb3WN3CVSoWvr2+GAvo3Q0NDqlWrxqlTpzRlgwcP5ty5c3h7e3Pq1Cnq1KnDo0ePNAthA3z66aeoVCo2bdqUYrvSMyhEyqKiojA3Nwfgww8/zDELor7Pk/AYas5QD85u/4ETP3xWWbcB5QRxUbC8ITy7oT4eHwIFDHUaksglbh5UTyhyb6TrSN4rP/QMPnnyhCpVqqCvr8/AgQMpXVo9zOPatWssWrSIxMREfH19k3WMpUWavwJMmzaN6dOnExwczIULF/Dz89P8ZGYiCODg4ICHh4dWWdmyZTX7Iad3uraRkZHWrCMLC4tMjVuI3MrMzEzz7/j06dOMGDGCq1ev6jiq9ytiaczcTyoBsl8xSYn/zBpOgirdwcgSyrWXRFC8X1IS/NwUNn4Gaz+CmDBdRyRQ3+08deoU5cuXZ8yYMbRv35727dszduxYypcvz8mTJzOUCEI6ksGXL19m28DyOnXqcP36da2yGzdu4OKivgWU0enaQojkSpYsiZmZGQBz587Fw8ODcuXK0adPHy5fvqzj6N7OxkK9lMSJm8/y55qDV3fB8kaw9UtIjAMjC6g9EAaeg49Xvf/5QujpQcvvIClBffxdCXUvodCpKVOmYGtry19//cWzZ884c+YMp0+f5tmzZ/z1119aE2vTK83J4CeffML+/fszfOHUGDZsGKdPn2bGjBncunWLDRs2sHz5cgYMUG+ZpFKpGDp0KNOmTePPP//k0qVLfPHFFzg6Or51oKUQ4t3Mzc3x8fFh0KBBFCyo3pHiypUrrFy5klGjRuk4urdzsjbRPP5o0d9ExOSjHUminsHv3eCRLwT8Ab9+pO4hBLCwB5lQI1LLsbJ63UkjK/WXivUdYcdAiHyq68jyrcmTJxMZGQlAwYIFqV69OjVq1ND8/5wZ0jxmcObMmcybN49WrVpRoUKFZNOdBw8enGnBAezatYsxY8Zw8+ZN3NzcGD58uGY2MaiXl5k0aRLLly/XTNdevHgxpUqVSvU18sOYAyHS49WrV2zevJk9e/bw22+/4ejoyMOHD3Ud1lstOXqb2XvVe4waFtDjwLD6uBQ203FUWSz8MSyrD1H/fFj33KteVFgWkRYZkRAHW3rD1T/Vx3bl4MsjUCD5Ys66lB8+v/X09AgODsbOzi7LrpHmZPBd3ZEqlSpDU5t1JT/8MgmREWFhYZp1tUqWLMmMGTP4+OOPdRvUW8zac42lx25rjo+NbJh3E8LA47CmzZvj6n2h1fe6i0fkLQlxcH4lHP8e6o+AGl/muC8Z+eHzW09PjydPnmit75zZ0rzGQGBgYFbEIYTIwaysrKhatSo+Pj7cvHmTTz75hBs3blCyZEldh5bM6BZlqO5akN5r1NtW+geF5d1k8NZBMDCD+CioOxwa/E/XEYm8pIAhfNgfqvcB/X/uAkY9A70CeXph6pyoVKlS710/9cWLF+luP90LTsXFxREYGIi7uzsFCsi6VULkdceOHWPNmjWaMbvXrl3LkckgQJOyRWhV0YHd/o/z5uzi20fAqQo0nQKNJ0D8KzBOvm+pEJlC/1/DwU4vhlMLoHIX9QLmZoV1F1c+MnnyZKysrLKs/TRncdHR0QwaNIg1a9YA6tm9xYsXZ9CgQTg5OTF69OhMD1IIoXtmZmZ8/fXX7Nu3jz///JM9e/bQpk2b9z9RRwqbqZdSOXztKT1qu+adXUliwmFdR+h3DIqUV39Q6799qyohMtVDH/XEEp/V6uWL2i7QdUT5QqdOnbJ0zGCaZxOPGTOGixcvcvToUa2VsD09Pd+6yLMQIu94PS5nyZIlfPrpp7x8+VLHEaXM1vzNUjPbL+TcSS9pEhkCs5xBSVTvN5xXElyRe3TeBJU6qx9f2gKxEbqNJx/Iji+yaU4Gt2/fzsKFC6lbt65WgOXKleP27dvveKYQIi/o06eP5vHmzZsZP368DqN5u/ZVnDSPh226yOFrT95RO5d4cPrNY9f6uotD5F8GxvDREihcUj1W9dIfuo4oz0vjPN90SXMyGBISkmJXZVRUVN65DSOEeKsPPviA8PBwze5AGzZsID4+563pV7SgKVv619Yc91p9ntshkTqMKIMUBW4eUD8u2QwajNRtPCL/Uqmgyhfqx0dnwquceXcgr0hKSsrSW8SQjmSwWrVq7N69W3P8OgH8+eefZdcPIfIJCwsLdu7cCUBoaCienp46jihlVV0KsqpHdc3xZ8tOk5SU9d+ys8TFjeCrHquNnce76wqR1Sp/DmZ2EPkEzv+i62jyNG9vb3bt2qVV9uuvv+Lm5oadnR1ffvklsbEZmyiX5mRwxowZjB07lv79+5OQkMD8+fNp1qwZq1atYvr06RkKRgiRe7i5udG5s3rs0PHjx+nSpUu23M5Iq0Zl7JjZoQKg3rf43N30L7+gU9f3vHlc+XPdxSEEqGcRd9kEnt9CnWG6jiZPmzx5stZWoJcuXaJ37954enoyevRodu7cycyZMzN0jTQng3Xr1uXChQskJCRQoUIF9u/fj52dHd7e3lStWjVDwQghcg+VSsWGDRvo3r07ABs3buTWrVs6jiplnWsUo4ZbIQA+W346dy4303QyNJ8F/U+Bbep3WBIiyzhVgbrD1HsaA9w5Bi/v6jSkvOjixYs0adJEc/zbb79Rs2ZNVqxYwfDhw/npp5/4/fffM3SNdC0Q6O7uzooVKzJ0YSFE3rB06VLNUlNBQUE5du3BvvWKczZQ3Sv429n7DGhUImePc37oCxfWq5ePqdYTChVXLwAsRE7k+yv8OUj9uExrqNoDHCqBedaOdcsPXr58SZEiRTTHx44do0WLFprj6tWr8+DBgwxdI809g3/99Rf79u1LVr5v3z727NmTwjOEEHmZsbExjRs3BtRLTIWFhek4opQ19ShCp+rOAHy//waLj97Okbe1AXUPy4pGcO5n2DUU9o3TdURCvJttWVD9s1XdtV2w/mP4viQcnaXbuPKAIkWKaHZ/i4uLw9fXlw8//FBzPiIiAgODjK01muZkcPTo0SQmJiYrVxRFFpwWIp+qWLEioJ71ZmNjw+HDh3UcUcr6N3RH75/OwO/2XWfx0Ry4HNaDc/Br2zfHXjOh1kDdxSNEajhXh0kvoPdBKOj6pvzoTO3xriLNWrZsyejRozlx4gRjxozB1NSUevXqac77+/vj7u6eoWukORm8efOmZkmJfytTpkyOHS8khMhas2bN0owdTEhIYPny5TqOKGUuhc04NfrN2JvjN0J0GM1/KAosrAEr/zUzu9c+qPU1WDroLi4h0sK5Ogy+AOOfQglPcKkLxRvqOqpcberUqRQoUIAGDRqwYsUKVqxYgaGhoeb8L7/8QrNmzTJ0jTQng1ZWVty5cydZ+a1btzAzy6ObwQsh3snIyIjVq1ezbNkyADZt2oS/v7+Oo0qZvZUxW/qrl8EKevlKd4H4roU9/1MngaBeuy3un3UQnapC951Q7MO3P1+InEqlggJG6sWpa/SFpOR3E0Xq2djYcPz4cV6+fMnLly9p37691vnNmzczadKkDF0jzRNI2rVrx9ChQ9m2bZumW/LWrVt88803tG3b9j3PFkLkZc2bN9c87t27N+fOndNhNG/nZG0KwOOwV0zZeYWqLgVpVTELe98UBQKPq/dzjY+GhFi4c0R9ztAMmkxUP246BZxrgHWxrItFiOxibgflPtJ1FHmGlZVViuWFChXKcNtp7hmcM2cOZmZmlClTBjc3N9zc3ChbtiyFCxfm+++/z3BAQojcq1ixYsyapR4w7uvrS1JSko4jSpmdhRGFzQxJUuCXvwMZsMGXOXuvkZgVC1IrChybrR4HeHkr3Nj7JhEEqNb7zeMKH0siKITIdiolHdPpFEXhwIEDXLx4ERMTEypWrEj9+rl3n8ygoCCcnZ158OABRYsW1XU4QuRq8fHxmJiYkJiYSJEiRejRowcTJ07E1NRU16Fpufc8ihM3nzF552XiE9X/Df7U+QPaVnLM/Ivd2AcbPgXHD6BSZzCyBH0DcG8Mphn/Vi9EfiWf35kjTcng6//kL1y4QPny5bMyrmwlv0xCZK7+/fuzdOlSzfGqVavo0aOH7gJ6h3vPo2jw3VEAKjhZsaV/bQwLpPmmibaEWLi2G8q0Uo+dCnuo3rbLqUrGAxZCaMjnd+ZI0/94BgYGFCtWLMWlZYQQ4rUlS5bw8OFDrK2tAZg7dy7x8fG6DeotXAqbMaejemmcSw/DmLXnWsYb9V4If/SE5Q3h7t9g5SSJoBAix0rzBJJx48YxduxY1q5dmymDFoUQeZOjoyNz586ld+/eBAQEMG/ePP73v//pOqwUeZW359fTdwl4GM4vfwdioK9iTMuy6Wss7CEcmqJ+HBEMdulsRwiR7w0fPjzVdefNm5fu66Q5GVy4cCG3bt3C0dERFxeXZMvJ+Pr6pjsYIUTe8umnnzJw4EBevXrFtWuZ0OOWRaxMDPhzQF3qf3eEoJev+NX7Hv9rXgY9vTRsV5cQp94t5ML6N2Wfb5YxgUKIdPPz89M69vX1JSEhgdKlSwNw48YN9PX1qVq1aoauk+Zk8KOPPsrQBYUQ+Ye5uTmLFy+mZ8+e/P777yxduhQjIyNdh5UiPT0VB4c3oMyEvbyKT2TFiTv0a5DKVf0jnsBPldXLxrxWf6R6vUAhhEinI0ferDwwb948LCwsWLNmDQULFgTU+xb37NlTa0eS9EjXbOK8RgagCpF1jhw5otm7uGvXrqxdu1bHEb2b1w/Huf4kAoBz4zyxtUhF8hr9Am7uV2+9VbwhNJkkPYJCZIP89Pnt5OTE/v37KVeunFZ5QEAAzZo149GjR+luO11T5kJDQ/n5558ZM2YML168ANRdlw8fPkx3IEKIvKlOnTpUrlwZgOPHj+s2mFRY0OUDzePAZ1Hvf0JcNBhbQ6VO6m242syXRFAIkenCw8MJCUm+hWZISAgREREZajvNyaC/vz+lSpVi9uzZfP/994SGhgKwdetWxowZk6FghBB5j6GhIX/99Reg/hafU2cVv1aqiAW1ihcGYPz2S8QmvGP1hBPzYHVLeH5Lvbi0Kg1jDIUQIg3at29Pz5492bp1K0FBQQQFBbFlyxZ69+5Nhw4dMtR2mpPB4cOH06NHD27evImxsbGmvGXLlrniW78QIvvZ29tjbGxMUlISDRs2zPEJYcki5gDceBJJ1akH+X7fdQ5dfcLTiBjtin5r4ZEfLKoO0c91EKkQIr9YunQpLVq0oEuXLri4uODi4kKXLl1o3rw5ixcvzlDbaU4Gz507R79+/ZKVOzk5ERwcnKFghBB5k0qlonXr1gCcOnUKQ0NDfvvttxy7Xd3wpqVo889OJJGxCSw8covea85TZ9ZhHrz4Z5LI1i/hxR314xE3wcxGR9EKIfIDU1NTFi9ezPPnz/Hz88PPz48XL16wePHiZCu7pFWak0EjIyPCw8OTld+4cQNbW9sMBSOEyLs2bdrEoEGDNMedO3dm2bJlOozo7axNDVnQ+QP+HFiHYZ6laOZRBID4RIV6c44wde5c8N+krlysFpjJ/31CiOxhZmZGxYoVqVixYoaTwNfSnAy2bduWKVOmaG7zqFQq7t+/z//+9z86duyYKUEJIfIePT095s+frxk/CPD111/zySefMHjwYH7++Wdy2uIGFYtaM8SzJMu/qMaKL6rhYGVMbb0AJkSoF5V+ZlQUevwlYwWFENnixIkTdO3alVq1amkm7a5du5aTJ09mqN00J4Nz584lMjISOzs7Xr16RYMGDShRogQWFhZMnz49Q8EIIfI2lUpFixYtuHnzJvr6+gD88ccfLFiwgL59+1KpUiU6dOjA559/zty5czM8Qy5TRL+A24dp6m6K95gmTPykLgB3kuzpH96L03df6jhAIUR+sGXLFry8vDAxMcHPz4/Y2FgAwsLCmDFjRobaTvc6gydPnsTf35/IyEiqVKmCp6dnhgLRpfy0TpEQOcXdu3c5evQowcHB/PLLL9y8eTNZnYIFC3L79m3NAqvZ6sE52DkEnl5WH3/1N9iXh9hIYh/4UmFVJHGJCiYG+lyY1BSjAvrZH6MQ+Vx++vz+4IMPGDZsGF988QUWFhZcvHiR4sWL4+fnR4sWLTI0byPNO5C8VrduXerWrZvuCwsh8jdXV1d69OgBwKhRo/jrr7949uwZMTEx7Nu3j+3bt/Py5UvOnTtHs2bNsje4syvgrxFvjk0Lg/LPZBcjc4xK1Ofn7iF88ctZXsUnsjcgmHaVnbI3RiFEvnL9+nXq16+frNzKykqzzF96pWvR6UOHDtG6dWvc3d1xd3endevWHDx4MEOBCCHyLz09PVq3bk2PHj346quv2LZtm2brSy8vL549e5a9Ad088OZxl99h5G1wqKhVpX4pW+qWUM8gHvLbBWb8dZVXce9Yk1AIITLA3t6eW7duJSs/efIkxYsXz1DbaU4GFy9eTPPmzbGwsGDIkCEMGTIES0tLWrZsyaJFizIUjBBCvFa7dm3N499//z17L37fW/1n161QyuutE0S+bvhm7+Llx+/wq/fdbAhOCJEf9e3blyFDhnDmzBlUKhWPHj1i/fr1jBgxgv79+2eo7TTfJp4xYwY//PADAwcO1JQNHjyYOnXqMGPGDAYMGJChgIQQAuCbb75hx44d/P3335w/fz57L+5aF67/BYXe/W27dgkbTo9pQvP5xwmNjsfnnkwmEUJkjdGjR5OUlESTJk2Ijo6mfv36GBkZMWLECK1lu9IjzT2DoaGhNG/ePFl5s2bNCAsLy1AwQgjxmp6eHp07dwZg1apVrFy5Mvsu7loPag2Egq7vrWpvZcysDupbyPuvPOHr9T6ERMRmcYBCiPxGpVIxbtw4Xrx4QUBAAKdPnyYkJISpU6dmuO10rTO4bdu2ZOU7duzQ7DAghBCZ4fW4QYBdu3Zl7cUUBSKfqh/X+hq8pqd6/cCaboWwMTcE4K9LwVSffpC5+6/nuHUThRC5V+PGjZk8eTKGhoZ4eHhQo0YNzM3NefnyJY0bN85Q22leWmbatGl8//331KlTh1q1agFw+vRp/v77b7755hssLS01dQcPHpyh4LJLfpqaLkRus3PnTtq2bQvA2bNnqV69euZfJCYc1nWEoLNQsRN0SPvOKIqiMG33VVaeDNSUbR9Qh8rO1pkYqBDi3/LT57eenh6FCxemTp06rF+/XrP7yJMnT3B0dCQxMf0T2NKcDLq5uaWuYZWKO3fupCuo7JaffpmEyG1evHiBvb29ZtejIUOG8MMPP6DKzF0//NbBjn/GO9f4Elp+l+6mgsNi+HDmIQDK2FuwfUAdjA1kDUIhskJ++vzW09PDz8+Pfv36ERUVxc6dO3F1dc2UZDDNt4kDAwNT9ZNbEkEhRM5WqFAh9uzZozmeP39+5v7/8uTym0TQvgI0m5ah5uytjBnYqAQA14IjmLzzitwuFkJkCgcHB44dO0aFChWoXr06R48ezZR207XOoBBCZKcmTZrw/PlzzfH333+fOQnW+V9gyZslbGgwGgoYZbjZvvWLa24Pbzx7nz5rzssahEKIDHl9N8TIyIgNGzYwZMgQmjdvzuLFizPctiSDQohcoVChQnTr1g2ApUuXcujQoYw3+tD3zeNWc6FMq4y3CViZGPBr7xoULWgCwKFrT6k8ZT8rjt8hMjaBuIQk6S0UQqTJf//PGD9+POvXr2fu3LkZbjvdexPnJflpzIEQudn9+/dxcXEB4Mcff2TIkCEZa/D2YXh2E5yqQtFqmRChtsQkhZGbL7LV72GK51tXdGBhlyqZfl0h8ov89Pl97949nJ2d0dPT7scLCAjAx8eH7t27p7tt6RkUQuQaxYoV4+uvvwbQum2cZuGPIDIE3BtDzX5ZkggC6OupmPdZZc6Oa8LARiUwN9Je53+X/2M+XnKKW08jeBoeI72FQoi3cnFxSZYIApQvXz5DiSCkYwcSIYTQJRsb9X7A6d6v+Ppe2PiZelHpch2gaNVMjC5ldhbGjPAqzbCmpYhPTCIuMYmuP5/BPyiM8/de4jnvOADGBnpUdy3EsKalqFKsYJbHJYTI2YYPH87UqVMxMzNj+PDh76w7b968dF8nXcngiRMnWLZsGbdv3+aPP/7AycmJtWvX4ubmRt26ddMdjBBCvM/rZPD48eMoipK2JWb8N8PWPurH3guheMPMD/Ad9PVU6OvpY2ygz5b+tVl5MpB1p+8RFh1PRGwCMfFJnLj5jBM3n9GotC0Lu1TBzEi+swuRX/n5+WmW1fLz83trvYwutZXm28RbtmzBy8sLExMT/Pz8iI1Vb7sUFhbGjBkzMhSMEEK8j62tLQCXL19m2bI0LA4dHwP7xr457n8KSjbN5OhSz0Bfj68auHPyf425NNmLa1Ob8/MX1figmDUAR66HsPrUXZ3FJ0RetWjRIlxdXTE2NqZmzZqcPXv2rXXj4+OZMmUK7u7uGBsbU6lSJfbu3atVJzExkQkTJuDm5oaJiQnu7u5MnTpVa9hHjx49UKlUWj8pbe37X0eOHMHa2lrz+G0/hw8fTt+b8Y80J4PTpk1j6dKlrFixAgMDA015nTp18PX1fcczhRAi45o2fZPA/f3336l/4oEJEPXPdnM9/oIi5TI5sowxNtDH06MIW/vXpkkZOwC+23edT5d6c+jqEx1HJ0TesGnTJoYPH86kSZPw9fWlUqVKeHl58fTp0xTrjx8/nmXLlrFgwQKuXLnCV199Rfv27bV66WbPns2SJUtYuHAhV69eZfbs2cyZM4cFCxZotdW8eXMeP36s+dm4cWOWvta0SPP9h+vXr1O/fv1k5VZWVoSGhmZGTEII8VaFCxfml19+oVevXmmbRHL3pPpPOw8o9mHWBJcJVCoVk9uV4/qTCIJevuLs3RecvfuCJmXsWNatKgX0Zd6fEOk1b948+vbtS8+ePQH1MlW7d+/ml19+YfTo0cnqr127lnHjxtGyZUsA+vfvz8GDB5k7dy7r1q0D4NSpU7Rr145WrdRLU7m6urJx48ZkPY5GRkbY29unKd73jRP872tLrzQng/b29ty6dQtXV1et8pMnT1K8ePF0ByKEEKlVsKB6csXLly9T94TEBPCaAWs/gp5/gV7O3h6uaEFTTv6vMRcfhNJ91VlCo+M5dO0p32y+yI+fVc7crfiEyCfi4uLw8fFhzJgxmjI9PT08PT3x9vZO8TmxsbEYGxtrlZmYmHDy5EnNce3atVm+fDk3btygVKlSXLx4kZMnTyZLzo4ePYqdnR0FCxakcePGTJs2jcKFC78z5neNE/y3jP6fkOZksG/fvgwZMoRffvkFlUrFo0eP8Pb2ZsSIEUyYMCFDwQghRGqkORnULwBu9aHucDDJPbN0Kzlb4zehKd1XneP4jRB2XHjEp9WcqVPCRtehCZGjREREEB4erjk2MjLCyEh7N6Fnz56RmJhIkSJFtMqLFCnCtWvXUmzXy8uLefPmUb9+fdzd3Tl06BBbt27V2gd49OjRhIeHU6ZMGfT19UlMTGT69Ol8/vnnmjrNmzenQ4cOuLm5cfv2bcaOHUuLFi3w9vZGX//tX06PHDmSpvchvdKcDI4ePZqkpCSaNGlCdHQ09evXx8jIiBEjRjBo0KCsiFEIIbQUKlQIUA9biY2NTfafPgBJibClt3riyAddoWxr8JyUzZFmnEql4vuPK1JjhnrHlc9/PsPeofUoY2+p48iEyDk8PDy0jidNmsS3336b4Xbnz59P3759KVOmDCqVCnd3d3r27Mkvv/yiqfP777+zfv16NmzYQLly5bhw4QJDhw7F0dFRs/5fp06dNPUrVKhAxYoVcXd35+jRozRp0iRNMV25coX79+8TFxenKVOpVLRp0ybdrzPNyaBKpWLcuHGMHDmSW7duERkZiYeHB+bm5ukOQggh0uJ1MgjQtWtXNm/enLyS7xq4vE392KGSOhnMpewsjfmp8wcM3qi+ZbTN7yFjWkgyKHKHo9ef4nPvJd80K51l17hy5QpOTk6a45S+INrY2KCvr8+TJ9oTsp48efLWsXy2trZs376dmJgYnj9/jqOjI6NHj9YaFjdy5EhGjx6tSfgqVKjAvXv3mDlz5lsXgy5evDg2NjbcunUr1cngnTt3aN++PZcuXUKlUmlmK7++Rfzv3sq0SvNI5F69ehEREYGhoSEeHh7UqFEDc3NzoqKi6NWrV7oDEUKI1HJ0dNQM1v7jjz/w9/dPXun57TePa/TNpsiyTttKjvSs4wrAw5evdBuMEKl0+NoTeqw6x4LDt7j6OPz9T0gnCwsLLC0tNT8pJYOGhoZUrVpVa1/zpKQkDh06RK1atd7ZvrGxMU5OTiQkJLBlyxbatWunORcdHZ1sZxB9fX2SkpLe2l5QUBDPnz/HwcEhtS+RIUOG4ObmxtOnTzE1NeXy5cscP36catWqcfTo0VS3k5I0J4Nr1qzh1avk/xG9evWKX3/9NUPBCCFEaqhUKnbu3KkZ+1OtWjXNwqwaUSHqP5tOAbO8Mcauuqu6R3SX/2N2XEh5v2MhcorbIZH0Wn0eAAcrY8wMdb+A+vDhw1mxYgVr1qzh6tWr9O/fn6ioKM3s4i+++EJrgsmZM2fYunUrd+7c4cSJEzRv3pykpCRGjRqlqdOmTRumT5/O7t27uXv3Ltu2bWPevHm0b98egMjI/7d33/FNVW0Ax39J2nQPRlso0JYhe5QtG1kFAZmKgLJcyJChIChDEEVQkamgL0sEQWQJMgRkU/ZeZa/SAYVu2qbJff+IXIgt2Ja26Xi+fvLpHSc3Tw61eXLuGbGMGDGCAwcOcP36dbZv306HDh0oU6YMAQEBaY49MDCQiRMnUrhwYbRaLVqtloYNGzJ58mQ++OCD56qXNP/LREdHoygKiqIQExNjMbrGaDSyceNGPD09nysYIYRIK41Gw6xZs3jttdcwGAzcvHmT0qVLm0+Gn4dTK8zbTnnn71Jlbzd1e8jyE9QrVQhPV/tnPEOI7PcwyUhsYjIDl5rnHrbRatg0pBHujnorRwbdunXj7t27jBs3jtDQUPz9/dm8ebP6xfLmzZsWrXwJCQmMGTOGq1ev4uzszMsvv8ySJUvUiaABZs2axdixYxkwYADh4eF4e3vz3nvvMW7cOMDcSnjq1CkWL15MZGQk3t7etGrVis8//zz1/s5PYTQacXFxAcy3vO/cuUO5cuXw9fUlKCjouepFo6RxZXStVvvMocsajYYJEybw6aefPldAz/LVV18xevRohgwZwvTp0wHzP9SHH37I8uXLSUxMJCAggO+//z7FaKFnuX37NiVKlODWrVsUL148i6IXQmQFPz8/bty4QWBgIC+++M/8gX9/Abunmrff3Qne1a0WX2Y7ExxFu1mPp7Wo7uPOsrdfxEGfs6fLEXlfgsHI9vPhjFp9ipiEZPX4vDdrElApffPrpVV++vxu1KgRH374IR07dqRHjx48ePCAMWPG8OOPP3L06FHOnDmT4WunuWVwx44dKIpCs2bNWLVqlUUHbr1ej6+vL97e3hkO5L8cPnyYefPmUbVqVYvjw4YN488//2TlypW4ubkxaNAgOnfunL6VCYQQuZaHhwc3btzg3r17jw9W7mJOBjvOzVOJIEDlYm6816QUf54K4faDhxy/Gcnxmw+oL9PNCCvaei6Md34+kuL4mLYVsiwRzG/GjBlDXFwcABMnTqRdu3Y0atSIQoUKsWLFiue6dpqTwSZNmgBw7do1SpQokaKzZFaKjY2lZ8+e/PTTT0yaNEk9HhUVxfz581m2bBnNmjUDYOHChVSoUIEDBw48biUQQuRZhQubk6C7d+8+PuheAgYcAM8KVooqa41uU4HRbSrQftZeTgdH8eWm8yx/tx7OdtbvkyXynwNXIywSwXcaleTVWiUoWdgJW1kxJ9M82b+wTJkyXLhwgfv371OgQIHnnnQ63f9Kvr6+aLVa4uPjuXDhAqdOnbJ4ZIWBAwfStm1bWrRoYXH86NGjGAwGi+Ply5fHx8fnqbOJCyHyFg8PD8A8IX7iN5Vg1Tugd8qzieCTmpQ1v/czwdFUHr+FYzfTOAm3EJnkx91XeP3HAwAUc3fg3MQAPm1bkbJeLpIIZoOCBQtmyopE6f4aeffuXfr27cumTZtSPf8889ykZvny5Rw7dozDhw+nOBcaGoper7foyAnm2cRDQ0Ofes3ExEQSExPV/ZiYmEyLVwiRvV544QXA/Lfn9/3X6OlR2soRZZ/3m5bmfEg02y+EAzB92yV+7lfHylGJ/OLPUyF8ufHxyh0/v1UHxxwwYjgvS0hI4NSpU4SHh6eYuuaVV17J8HXT/a82dOhQIiMjOXjwIE2bNmXNmjWEhYUxadIkvv322wwHkppbt24xZMgQtm7dmmJtwOcxefJkJkyYkGnXE0JYz4cffsj3c2YTGhZOUIQRXPJP/yQnOxt+6lWLH3Zd4estQey+eJfRq0/xTqNSlPKQhQBE1rkZEc/w304AUNrDiT8/aIS9rQxiykqbN2+mV69elv2j/6HRaLJ30um///6badOmUatWLbRaLb6+vrzxxhtMnTqVyZMnZziQ1Bw9epTw8HBq1KiBjY0NNjY27Nq1i5kzZ2JjY4OXlxdJSUlERkZaPO9Zs4kDjB49mqioKPVx7ty5TI1bCJF9HB0dGdrZPGHs57uTOP4w7ZO45gVarYZ+DUri7mgLwK+HbtHs210EfLebvZfukcYJI4RIl18P3yQx2UQ5Lxc2D20siWA2GDx4MK+++iohISGYTCaLx/PelU13MhgXF6fOJ1igQAG103aVKlU4duzYcwXzb82bN+f06dOcOHFCfdSqVYuePXuq27a2thaziQcFBXHz5s1nziZuZ2dnMVP5o3l7hBC5UMhJ6j/cqu4uPJj/+s056HVsGdqY/k1KU7yAAwBBYTG8Mf8gH686RbLx6SshCJFeBqOJlUduAzCs5QvSNzCbhIWFMXz48HRNnZdW6b5NXK5cOYKCgvDz86NatWrMmzcPPz8/5s6dm65lVdLCxcWFypUrWxxzcnKiUKFC6vG33nqL4cOHU7BgQVxdXRk8eDD16tWTkcRC5BdHF9HI14YP6+n5NjDpmf2F8zIvV3tGtSnPyIBybDwTwqBl5nWMfztym/txSfyvd20rRyjygmv34nj/l6Pci02ksLOeZuUzPzERqevatSs7d+58PLl+Jkp3MjhkyBBCQkIAGD9+PK1bt2bp0qXo9XoWLVqU2fH9p++++w6tVkuXLl0sJp0WQuQTLw4Ar0r4252CwFk8eJD/WgafpNVqaFfVm5YVvXht3gFO3opk98V7GE0KOu3zjzoU+ZPBaGLKpgv8b+819Vj3Oj7obaRVMLvMnj2bV199lT179lClShVsbW0tzj/PknRpXoHkaR5NMePj46PO95Xb5KcZzIXIq/7880/atWtHjRo1OHr0qLXDyRGMJoUK4zaTlGyianE31gxoIAmhSDeD0US3eYEcuxkJQGFnO4a3LEvXmsWtngzmp8/v+fPn079/f+zt7SlUqJDFlDIajYarV69m+Nrp/lecOHEi8fHx6r6joyM1atTAycmJiRMnZjgQIYRIt3N/wL4ZgLkPM8CxY8eYNm2aNaPKMXRaDS+VM89FeOp2FIev37dyRCI3mr/3mpoIvl67BAc/aU6PutIqmN0+/fRTJkyYQFRUFNevX+fatWvq43kSQchAy6BOpyMkJEQdRPJIREQEnp6emT7PYHbIT98shMgzws7BD/WgdDPo8RtRsfEWc45GRERYLJuZX5lMCo2m7iA48iEAtf0K8GWnKrzgJQPnxLPdjUnkiz/PsfbEHQC+ebUaXWvmrM/I/PT5XbBgQQ4fPpwlfQbTndYripLqbNcnT56UP7xCiOyRFAfzGpm375wAjQ43Nzfu3LmjFqlUqZLF5PL5lVar4dO2j1djOXz9AT/teb5WBJH3XQ6PpfeCQ2oi2KKCJ11qFLNyVPlb7969n3sN4qdJ8wCSR2vfaTQaypYta5EQGo1GYmNj6d+/f5YEKYQQFiJvgSnZvN3lJ/hnrfSiRYsyduxYPv/8c0JDQ/H29ubKlSspVinKb16uUpRLX7Rh3Loz/HroFr8duY2iwKu1SlDKw4nCznbWDlHkIGPWnuaXAzcBcHOw5ePW5XmtVvFMWfZMZJzRaGTq1Kls2bKFqlWrphhA8jzdY9KcDE6fPh1FUejXrx8TJkzAzc1NPafX6/Hz83vm3H5CCJFpYsPMPwuXgzKWa5ZPnDgRFxcXRo4cyf3799m7dy/t2rWzQpA5i61Oy/CW5fj96G0MRoWVR2+z8uhtbHUa/hrWhJKFnawdosgBlh+6qSaC5Yu48MMbNeV3I4c4ffo01atXB+DMmTMW5543UU9zMti7d28ASpYsSYMGDbCxkfUHhRBW8igZdPZM9fSIESNYtWoVBw8eJDw8PBsDy9k8XOw4/VkA3++8wsGrERy8dh+DUWHQsmMs6lsHDxdpIczPYhIMfLrWnGSULOzEpiGNpDUwB9mxY0eWXTvdfQZdXFw4f/68ur9u3To6duzIJ598QlJSUqYGJ4QQqbr0l/mny9Mnun/Uyfrfy1Xmd/a2Ooa3LMuK9+oxIqAcAGfvRPNz4HXrBias7uiNBxhN5jGlS96qI4lgDmIwGGjevDmXLl3KkuunOxl87733uHjxIgBXr16lW7duODo6snLlSkaOHJnpAQohRAqNPgStDZRv+9Qij/oJ5vdJqJ/l/SalqeHjDsDPgTcIjUqwbkDCKk7djuStRYfps/AwAM3Ke1K8gKOVoxJPsrW15dSpU1l2/XQngxcvXsTf3x+AlStX0qRJE5YtW8aiRYtYtWpVZscnhBApeZSHUbegUsenFnk07+CkSZMICgrKpsByF61Ww/CW5tbBqIcGXv8xkARD7pseTGScoiiMWnWa7RfM3Sm83ex5u1FJK0clUvPGG28wf/78LLl2ujv+KYqCyWRe9Hzbtm1qx+wSJUpw7969zI1OCCGeZDLBrYNQog7on91yUadOHXV7yZIlTJo0Kaujy5XqlS7EWw1LMn/vNa5HxFN+7GZerVmczztWxt5WZ+3wRBY7cuMB50KiAfiqcxW61iyOjU4mk86JkpOTWbBgAdu2baNmzZo4OVkO7Hme0cTp/hevVasWkyZNYsmSJezatYu2bc23aa5du4aXlyxYLYTIQofmwdKusPodCDn5zKKvvPKKulbnF198weLFi7MjwlxHp9Uwtl1F+jd5PJHtyqO3eWX2XrltnMclGIxM+tM8BqBbrRK8XsdHEsEc7MyZM9SoUQMXFxcuXrzI8ePH1ceJEyee69rpXoHk1KlT9OzZk5s3bzJ8+HDGjx8PwODBg4mIiGDZsmXPFZA15KcZzIXIlRQFTq2ANe+Z9919of8esHd75tPOnTtHpUqVAPDx8eHq1avodNLa9TThMQkM+fUEgVcj1GMtK3rxRafKeLrYWzEykZkURWH7+XBm77jMiVuR2Ntq2TC4IWU8c9+qNPL5nTnSnQw+TUJCAjqdLsUkiLmB/DIJkYMlJ8KPL0H4WfO+1gaGnAK3tK2GsGfPHho3bgxA06ZNs3R6hrzAYDQxbetFfth5RT3mW8iRnR81ldGlecS6E8EMWX4CAGc7G/7XuxYvlipk3aAyKL99fkdGRjJ//nx1VpdKlSrRr18/i7mfMyLT2oPt7e1zZSIohMjhji56nAh614C+m9KcCAI0aNCAbt26AbBz504ZXfwfbHVaPm5dnjMTAuhR1weAGxHx9FpwSJ12RORum8+EAlDdx501A+rn2kQwvzly5AilS5fmu+++4/79+9y/f59p06ZRunRpjh079lzXls4BQoicLfyfeU3LtoZ3d5gHj6SDVqtl+fLllChRAoD169dndoR5krOdDV92qkKTsh4A7Ll0j2M3JZHO7UwmhYPX7gMwpm0FXvDKfbeG86thw4bxyiuvcP36dVavXs3q1au5du0a7dq1Y+jQoc91bUkGhRA5W9tvoVJnaDj8uS5TuXJlwLyaUmhoaGZEli/M7F5d3T5yXZLB3O58aDT345JwsNVRpZi7tcMR6XDkyBE+/vhjixXgbGxsGDlyJEeOHHmua0syKITI2bQ6eHUh+NR9rss8OSn+oUOHnjeqfMPNwZbXa5tbVadsvsC+yzKFWG517OYD3l5sThpq+RVAbyMpQG7i6urKzZs3Uxy/desWLi7P18Kb7t+EiRMnEh8fn+L4w4cPmThx4nMFI4TIpxQFru2Bc3/AmVWwYzJ85gYr3nx8m/g5NW3alFdffRWAd999N9U/qiJ1nao/7qO56uhtK0YiMiLBYOTVufvp/P1+QqIS0GqwmEpI5A7dunXjrbfeYsWKFdy6dYtbt26xfPly3n77bbp37/5c1073aGKdTkdISAienpYLxEdERODp6YnRmPtmr89vo5GEyFEirsD8lhAfkfKcRgsfHIcCfpnyUgsXLqRfv34AVKtWjd27d+Pq6pop187rlh68wadrzgAwoGlpRrYub+WIRFr9sPMKUzZfAKB+6UJ88nIFKhd7vtGnOUV++vxOSkpixIgRzJ07l+TkZMC8TN3777/PV199hZ2dXYavne6WQUVRUp1e4OTJkxQsWDDDgQgh8on4+3BwHsQ80W/vUSLoWQn8GkHp5tBgCPT7K9MSQYA+ffowbNgwwPw3a+zYsZl27byuWfnHDQDf77zCon3XrBiNSKtdF++qieD49hVZ9s6LeSYRzG/0ej0zZszgwYMH6kTT9+/f57vvvnuuRBDSsRxdgQIF0Gg0aDQaypYta5EQGo1GYmNj6d+//3MFI4TIw+LuwZZPzJNHA3hVBpci4FIUXvvZvF8oa29daTQaJk2axNq1a7l27dpzz9qfnxR1c+DSF22oMXErMYnJfLb+HK/X8ZEl63K473dcBsBRr6NXPT/rBiMyhaOjI1WrVs3Ua6Y5GZw+fTqKotCvXz8mTJhgMcGhXq/Hz8+PevXqZWpwQog8ZEknCD31eP/RF0q9I1TskG1hODo6snDhQpo2bUpISEi2vW5eYKvTsmloIxpOMU/cXX7sZioWdcXZ3obEZBO1fQswtGVZnO3Svey9yALX7sWp08j8NawxOq1MGp4babXa/5zwXaPRqLeOMyLN/8f27t0bgJIlS9KgQQOLoc1CCPFMiTGPE8EGQ6DZONBZ729I0aJFASQZzIDiBRzp36Q0c3eZVyg5FxKtnjt5K5L/7b1GaQ8nnOxsGNW6PPXLFLZWqPmOoiiExyQSFp3AqdtRzNx+CYCXynlQvICjlaMTGbVmzZqnngsMDGTmzJmYTKbneo10DyDZuHEjOp2OgIAAi+NbtmzBZDLRpk2b5wrIGvJTB1Qhst3DB/BtBUh+CC7e8GHmjA5+HjExMerAkSpVqtCsWTO+/fZbWbc4HaITDOwKuotJUdBpNfwceIND/7RCPalLjeKMa1cRN0dZoSorXAyLYcyaM1wMj8HeRkdodEKKMgv71Oal8p6pPDv3y6+f30FBQYwaNYr169fTs2dPJk6ciK+vb4avl+6v5qNGjeKrr75KcVxRFEaNGpUrk0EhRBa6d9ncN/DBNfDvYe1oAHBxcaFZs2b8/fffnD59mtOnT+Pg4MAXX3yBVitzr6WFq70t7at5q/vtqnrzIC6JfVfuYTQp6tq3q47dZtWx25TycCLRYCI+KZk3X/RleKtyVoo8b1AUhenbLjHjn9Y/MwMABZ30lCjoSNVibrxepwSVvGXASF5x584dxo8fz+LFiwkICODEiRPqhPrPI90tgw4ODpw/fx4/Pz+L49evX6dSpUrExcU9d1DZLb9+sxAi25hMkBgFDgWsHYlKURRu3rzJoEGD2LBhA2C+HdOxY0frBpZHJBiMDFp2nG3nw1I93+iFwix56/kmEs+PFEXhQmgMKw7fYtH+6+rxwc3K0KSsBxW9XXHU559uXPnl8zsqKoovv/ySWbNm4e/vz5QpU2jUqFGmXT/dvzFubm5cvXo1RTJ4+fJlnJycMisuIURulxQHt49Aibpga5+jEkEwd7j29fXlm2++UZPBw4cPSzKYSextdfyvdy0SDEYOX7+PjVaLva2WoStOcCMinj2X7rHpdAhtqhS1dqi5xsJ915ix/RKR8Qb1WNNyHvz4Zi1ZTSQPmzp1KlOmTKFIkSL8+uuvdOiQ+QPu0t0y+N577xEYGMiaNWsoXdo8DcTly5fp0qULtWvX5n//+1+mB5nV8ss3CyGyzebRcOB70DtD7begZh8oWMraUT3Vt99+y0cffQTA6dOnM+W2i0idyaTQaOoOgiMfAvDbe/WoU9JyjtpTtyPZcjaU9tW8KV9EJgUH88jgFtN2YTQp2Gg1VPJ2pVc9P7rUzN+fWfnh81ur1eLg4ECLFi2e2a959erVGX6NdLcMTp06ldatW1O+fHm14m/fvk2jRo345ptvMhyIECKPCD5mTgQBdHq4exHscnafperVq6vbH374IVu2bLFiNHmbVqvh61er0uOngwC8Ni+Qcl4udKlZjDaVi7Jg3zUW7rsOwJwdV/Ar5MicnjXyfb+3GdsuYjQp1PYrwKK+dXCS6XvyjV69ev3n1DLPK90tg2Dus7B161ZOnjyJg4MDVatWpXHjxlkRX7bID98shMgWcRHw9T8tgM5eMPw8aHP+CF1FUZg9ezYffPABAAMGDGDOnDlWjipvC7wSQc//HcCUxk8gv0KOdK/jg09BR2r6FsDT1T5rA8xBgkJjCJi+G4ANgxvKCiJPkM/vzJGhZPCRhIQE7OzssjxjzWryyyREJoi9CwvbQMQ/oxtf+zlbJ5N+XoqiUKVKFc6ePYuLiwtRUVG5/m9bTpdsNHEuJJqhy09w9Z558GH5Ii5ULe7GR63K8cfJO0z6M/WpiP78oGG+aS2c/fclvvnrIs3Ke7KgT21rh5OjyOd35kh3O7PJZOKLL75g7ty5hIWFcfHiRUqVKsXYsWPx8/Pjrbfeyoo4hRA5nd4Rqr4G+2ZCjTdzVSII5gElhw8fxtHRkZiYGO7fv0+hQoWsHVaeZqPTUrW4O39/1JQEg5GIuCSKuTuo599uVIrudXxYcfgWm8+EEhGXyJW75qSx7cy9dKtVggpFXbgfl0SlYm40K++JrS7vDaTYc+keYJ48WoiskO5kcNKkSSxevJipU6fyzjvvqMcrV67M9OnTJRkUIj9SFNA7Qf0PoHBZKNPc2hFliIODA4ULF+bevXtUqlSJ4OBgmYg6m9jb6iwSwUec7Gzo17Ak/RqWBGDx/uuM/+MsACuO3LIoW8bTmS1D886yazEJBlYcvqUuKdfwBUkGRdZI91eon3/+mR9//JGePXta/JGsVq0aFy5cyNTghBA5nKLAvhmwYRgkxZunkKnUEexcrB1ZhvXoYZ4YOywsjMuXL1s5GvFvvev7sWfkS/Sp70eLCp40f2Jljcvhsby9+DBJyc+3NJc1KYrCgasRbD0XRsB3u9Xb5CUKOuBXSJaUE1kj3S2DwcHBlClTJsVxk8mEwWBI5RlCiDzJmAwLAiD4iHn/4mYYchJs7Kwb13OaMWMGf/75J1euXCE4OJhy5WSljJymREFHPnulksWxN+cfZM+le+wIukuVz7bwx6CGlCuS+76ULDt0k0/XnFH3HfU6Ovh706VGcenDKrJMupPBihUrsmfPnhRr4P3+++8W0zMIIfIwYzJsHfs4ESz1EvT8HXR5Y7oLHx8frly5wh9//EGzZs2sHY5Igy86VqH7TwcIjnxIYrKJDnP2svTtulQs6oaDPmtu9T9MMvLVpvOsOHILL1d74pOM2Gg1FHGzp1RhZ14sVZAXSxWieAEHQqIScHWw5WxwFHFJyfgUdKSMp2Wyejk8hs83nFP3Czja8l03f5qWy5vrCoucI91/uceNG0fv3r0JDg7GZDKxevVqgoKC+Pnnn9VZ/IUQedi13fBLFzAmmfdLN4c3Mz7ZaU70aFTijBkzAPP8qnq93pohif/gU8iRXSOa8vVfQczbdZUEg4kuPwQCsKBPLZqV9wLMCVfglQiu3I3jTHAUd2MT8XKxp2P1YnSuUQx722cnjgajiVv34wmNTmDB3mtsOx8OwI2IeLVMSFQCx29GsurY7TTFXszdgRq+BVh/8g4ADcsU5ud+dQDzvIxCZLUMTS2zZ88eJk6cyMmTJ4mNjaVGjRqMGzeOVq1aZUWMWU6GpguRDmsHwIml5u1qPaDZp+CWt/6/OXnyJA0bNiQ2NhYAGxsbmjVrRocOHRgwYICVoxPPoigKc3ZcZtfFuxy+/kA93sHfG51Gw7qTdzA+Y3LDPvX9aFnRi9p+BVMs8RYVb6DT9/vUaXCe9FbDklwMi+Hg1fskGTPeZ1Gv07L9wyaUKCj9A9NCPr8zR5qSwZkzZ/Luu+9ib2/PzZs3KVGiRJ7quyC/TEKkQ+xduHMMHAtD8ZrWjibL3L9/n3feeYe1a9diMj3+cF+zZo2sX5xL7Lt8j57/O5jieL1ShShf1IWShZ24GBbDLwdupvr8Ri8UZkDTMtT0LcD+K/eYv/eaOs0LgK1Ow6zuNWhduUiK5yqKwkODka3nwngQl8SLpQtRzN0BvY2Wy+GxXA6P5a+zYfx1LhSD0fwxrNNq+LJTZbrV9smkGsj75PM7c6QpGbSxseHOnTt4enqi0+kICQnB0zPv9GGQXyYh/oPJBEcXgE998KwAeejL4H+Ji4sjMDCQ1q1bYzQaAfDz82P+/PnSnzAXuHU/nkl/mvvhuTvoaVOlSKp98E7djmTcurOcuBX5zOvZaDX81KsW50KiaVXRixe8ct8glbxEPr8zR5qSQR8fH0aPHs3LL79MyZIlOXLkCIULF35q2dxGfpmESIWiwO3DcOFPuHMcru0CvQt8fD3PDBRJj0OHDjF8+HD27dunHvvss88YP368FaMSmS06wcDWs2F89sdZYhKTLc55u9kzu2cNavgUsFJ04t/k8ztzpCkZ/PHHHxk8eDDJyclPLaMoChqNRv3mnJvIL5MQT1AUWDcIzq2FpFjLc+Vehu6/WiWsnOLkyZP4+/sD5oEmN2/ezFPdZsRjUfEGdl26S2R8EnVKFqSMhzM2eXCFk9xMPr8zR5oHkMTExHDjxg2qVq3Ktm3bnrpMU7Vq1TI1wOwgv0xCPOHuRZjzxPqnhcuBXwNw94Wq3cC1qPViyyHCwsIoUsTcT6xNmzZs3LjRyhEJkT/J53fmSPNXHBcXFypXrszChQtp0KAB1apVS/UhhMjlCr8AH5yAgMkw6hYMOgTtvoOGQyUR/IenpyedOnUCYNOmTdy9e9fKEQkh0mrOnDn4+flhb29P3bp1OXTo0FPLGgwGJk6cSOnSpbG3t6datWps3rzZoozRaGTs2LGULFkSBwcHSpcuzeeff86TbW2KojBu3DiKFi2Kg4MDLVq04NKlS1n2HtMr3e3dvXv35uHDh/zvf/9j9OjR3L9vXjPx2LFjBAcHZ3qAQohs9OAGJEZDwZJQbwDYu1o7ohxJo9GwevVqtXVwzZo1Vo5ICJEWK1asYPjw4YwfP55jx45RrVo1AgICCA8PT7X8mDFjmDdvHrNmzeLcuXP079+fTp06cfz4cbXMlClT+OGHH5g9ezbnz59nypQpTJ06lVmzZqllpk6dysyZM5k7dy4HDx7EycmJgIAAEhISsvw9p4mSTidPnlQ8PDyUMmXKKDY2NsqVK1cURVGUTz/9VHnzzTfTe7kc4datWwqg3Lp1y9qhCGEdZ9cpyrTKijK1jKKcWaMokbetHVGu0KZNGwVQAOX7779XjEajtUMSIl9J7+d3nTp1lIEDB6r7RqNR8fb2ViZPnpxq+aJFiyqzZ8+2ONa5c2elZ8+e6n7btm2Vfv36PbWMyWRSihQponz99dfq+cjISMXOzk759ddf0xR3Vkt3y+CwYcPo06cPly5dwt7eXj3+8ssvs3v37kxLUoUQ2WT3N/DbmxB1E+LC4e9JEH/vv58nGDVqlLo9YMAAFixYYMVohBDPkpSUxNGjR2nRooV6TKvV0qJFCwIDA1N9TmJiokWuA+Dg4MDevXvV/fr167N9+3YuXrwImAeZ7d27lzZt2gBw7do1QkNDLV7Xzc2NunXrPvV1s1u6k8EjR47w3nvvpTherFgxQkNDMyUoIUQ2MSbDzsnmbXs3GH4eBh+BotL/Ny0aN27M6dOn1f2DB1NOcCyEyHoxMTFER0erj8TExBRl7t27h9FoxMvLy+K4l5fXU/OXgIAApk2bxqVLlzCZTGzdupXVq1cTEhKilhk1ahSvv/465cuXx9bWlurVqzN06FB69uwJoF47Pa+b3dKdDNrZ2REdHZ3i+MWLF/Hw8MiUoIQQ2SQuHEz/TBk15BS4els3nlyocuXKaovg0qVLiYqKsnJEQuQ/FStWxM3NTX1Mnjw5U647Y8YMXnjhBcqXL49er2fQoEH07dsXrfZx+vTbb7+xdOlSli1bxrFjx1i8eDHffPMNixcvzpQYskO6k8FXXnmFiRMnYjAYAHNH6ps3b/Lxxx/TpUuXTA9QCJGFYkJBpweXouDgbu1ocq2SJUsC8PDhQ0qXLs2RI0esHJEQ+cu5c+eIiopSH6NHj05RpnDhwuh0OsLCwiyOPzlV1L95eHiwdu1a4uLiuHHjBhcuXMDZ2ZlSpUqpZUaMGKG2DlapUoU333yTYcOGqQnpo2un53WzW7qTwW+//ZbY2Fg8PT15+PAhTZo0oUyZMri4uPDFF19kRYxCiMxmeGheYs6zItTqZ55ORmRYw4YN6devHwARERHUrl2bAQMGWEwtIYTIOi4uLri6uqoPOzu7FGX0ej01a9Zk+/bt6jGTycT27dupV6/eM69vb29PsWLFSE5OZtWqVXTo0EE9Fx8fb9FSCKDT6dQ1zUuWLEmRIkUsXjc6OpqDBw/+5+tml3SvKeXm5sbWrVvZu3cvp06dIjY2lho1alh0jBRC5GD7Z8FfY2DoGXAvAdW6g1Pqy0uKtLGxsWH+/Pk0bNiQd999l+TkZH744Qd69+5N3bp1rR2eEOIfw4cPp3fv3tSqVYs6deowffp04uLi6Nu3LwC9evWiWLFiaqvewYMHCQ4Oxt/fn+DgYD777DNMJhMjR45Ur9m+fXu++OILfHx8qFSpEsePH2fatGnqF0SNRsPQoUOZNGkSL7zwAiVLlmTs2LF4e3vTsWPHbK+D1GR4gdGGDRvSsGHDzIxFCJGVFAXOrDInggD3r5qTQW9/q4aVl/Tt25c2bdpQtKh5cu7AwEBJBoXIQbp168bdu3cZN24coaGh+Pv7s3nzZnVwx82bNy1a+RISEhgzZgxXr17F2dmZl19+mSVLluDu7q6WmTVrFmPHjmXAgAGEh4fj7e3Ne++9x7hx49QyI0eOJC4ujnfffZfIyEgaNmzI5s2bU4xUtpY0L0cH5ubURYsWsXr1aq5fv45Go6FkyZJ07dqVN998M9euzynL2Yg8LykO5jWBiCdmvP8wCFxyRn+VvGbYsGFMnz4dgJ07d9KkSRPrBiREHiWf35kjzX0GFUXhlVde4e233yY4OJgqVapQqVIlbty4QZ8+fdSlmYQQOYyimOcSfJQIlngRev0hiWAW6tGjh7rdqVOnp65uIIQQOUGabxMvWrSI3bt3s337dl566SWLc3///TcdO3bk559/plevXpkepBDiOZxbC3unmbcrdoTXcs90B7lV7dq12bRpE23atOHBgwfUrl2bI0eOyPRbQogcKc0tg7/++iuffPJJikQQoFmzZowaNYqlS5dmanCTJ0+mdu3auLi44OnpSceOHQkKCrIok5CQwMCBAylUqBDOzs506dIlxfBtIfIVowF+7gDTKpr3y7SEXuvArxE0Gfns54pM06pVK/XL8c2bNylbtiwrV66UEcZCiBwnzcngqVOnaN269VPPt2nThpMnT2ZKUI/s2rWLgQMHcuDAAbZu3YrBYKBVq1bExcWpZYYNG8b69etZuXIlu3bt4s6dO3Tu3DlT4xAi10iIhrmN4OpOiA42H7NzhlJNoc8G8KpkzejyFa1Wy6JFi9RRh5GRkbz22mt4enpSq1YtGjZsqPYrFEIIa0rzABK9Xs+NGzfUUXL/dufOHUqWLJnqEjCZ5e7du3h6erJr1y4aN25MVFQUHh4eLFu2jK5duwJw4cIFKlSoQGBgIC+++GKarisdUEWekBQPv74O13aZ98u9DN1/tW5MAoCzZ88yfPhw/vrrrxTnfv31V15//XUrRCVE1jOZTEyaNIkDBw7w559/ZvpAU/n8zhxp7jNoNBqxsXl6cZ1OR3JycqYE9TSPlnkqWLAgAEePHsVgMFjMcVi+fHl8fHzSlQwKkSfsn/k4EazYEbousGo44rFKlSqxZcsWbty4wc2bN4mJiaFt27YAdO/enaVLl+Lj40OFChVo164dfn5+1g1YiEzy22+/MX78eADWr1/PK6+8YuWIRGrSnAwqikKfPn1SndUbyNIWQTB/uxg6dCgNGjSgcuXKgHnxZ71ebzHfD/z34s+JiYkW8cbExGRJzEJkq73fmX+6+cArs0Crs248IgVfX198fX0BOHDgAI0aNcJgMLBhwwa1zODBg2nVqhUbN25Ep5N/Q5F7JSYm0r17dwBq1qxJmzZtrByReJo0J4O9e/f+zzJZOZJ44MCBnDlzhr179z73tSZPnsyECRMyISohrMTwEE4uh4odwNHcUk6ZFnDhT3hjFdi7Wjc+8Z/q1q1LeHg4y5cvJyYmhhs3bjB37lyMRiN//fUXrq6uXLx4kWLFilk7VCEyZMSIEer2unXrsLW1tWI04lnSNem0tQwaNIh169axe/dudUF4ME9p07x5cx48eGDROujr68vQoUMZNmxYqtf7d8tgcHAwFStWlD4HIme6dRj2TQeHAuDkAZe3Quhp87lyL8Mrs8GpEFzbDTo9+Ej3iNzKaDQyYcIEPv/8cwCmTp1q8YEqRG6wbds2tm3bxpQpUwDo3Lkzq1atypLXkj6DmSPDy9FlB0VRGDx4MGvWrGHnzp0WiSCYm51tbW3Zvn07Xbp0ASAoKIibN28+c/FnOzs7i9vd0dHRWfMGhMgMEZfhwobUzxniIT7CnAyWbJy9cYlMp9PpmDhxIhqNhokTJzJy5EhOnjzJyJEjqVq1qrXDE+KZTCYTCxcu5O2331aPVaxYkd9//92KUYm0yNEtgwMGDGDZsmWsW7eOcuXKqcfd3NxwcHAA4P3332fjxo0sWrQIV1dXBg8eDMD+/fvT/DryzUJkiTOrwLMieJSH1EbQJUTD4Z/MrXmOhUDvZH7YOED4OSjqDyVqw73LcG0nxD+A+HsQf9+8pnCjD83lRZ5z+fJlKlasiMFgUI/99NNPFh+yQuQkhw8fpk6dOhbH/P39mTZtWqrzE2cW+fzOHDk6GXzaEPSFCxfSp08fwDzp9Icffsivv/5KYmIiAQEBfP/99xQpkvaltuSXSWS6pDiYWgqSE8CnPjQZAaWbmc8pijk5TIyBzaPg+C+pX8OzIvj3gPqDsy9ukWMoisL8+fN555131GNXrlyhVKlSVoxKiJQURcHf359Tp04BMHbsWMaPH58tA6Dk8ztz5OhkMLvIL5PIdJE34c+P4PI2UIyPjxd6Adx94M3V5v34+7DtM4i6bb7lmxgLyQ/NCaO3P9QbBMVqWOMdiBzi8uXLvPDCC+p+eHi4LGsncgxFUXj99df57bffADh//jzly5fPtteXz+/MkaP7DAqRq9y7bJ70ucqrUKUr9PwNws7C7m/g7BpAgYhL5kfwMXOS51gQXplp7chFDlamTBnmzZvHe++9B8DGjRvTNLuDENmhS5curFmzBoBPP/00WxNBkXmkZRD5ZiGeQ/AxuHUQ4u7Bnm8eH++zEfwaPN6/dRjO/wFFqkDx2uDuC9o0rwYpBH369GHx4sUA7Ny5kyZNmlg5IpHfrV27lk6dOgEwceJExo4dm+0xyOd35pCWQSEyKiYMfkqlY3Sdd8HnX6PZS9Q2P4TIoN69e6vJ4Pr16yUZFFZ1//59dUBTmzZtrJIIiswjTRNCZJTNP9MT6ezAuYj5Z5OP4eWvpdVPZLqXXnqJr7/+GjCvBS+ENa1bt46IiAiKFCnCypUrrR2OeE7SMihEet05AQVLmad1GXwMCpW2dkQin/D29gbg119/ZcKECRYDS4TITmvXrgXM07s5OckUV7mdNF8I8SxRwXDv0uP90NPwYxP4qoR59K8kgiIb+fn5qdtDhgyxXiAiX4uLi+Ovv/4CoGPHjtYNRmQKSQaFeJLRAFf+hrgI837ISZhdC2bWgLmNYG5D83GXomDvZr04Rb5Ur149dd7BTZs2MXHiRC5duvQfzxIi80RHR9OwYUMSEhIoWbIkVapUsXZIIhNIMijEk/ZMgyWd4O4F875iBDRw/wqEnnpcroG0yojsp9FomDlzJi4uLgCMHz+eihUrEhwcbOXIRH4QGBhIkSJFOHHiBABdu3Z96uIQIneRPoNCPHL6d9j5pXn79iHz1DAV2sPH1+D6Pnj4AGLDoHJnc59BIazA3t6e3bt3s379esaNG0dycjInT56kWLFi1g5N5GG//PILb775prrfr18/xo0bZ8WIRGaSZFCIR3ZNfbxducvjbYcCUKFd9scjxFP4+/vj7+/P4cOHWb9+PRMmTKBJkybSkV9kiZCQEItlEZcvX063bt2sGJHIbHKbWIhHov+51dZtqXnJOCFyuEf9tQ4dOkTx4sW5cuWKlSMSedFXX31FQkICJUqU4ObNm5II5kGSDAphNMDKvpAUa94v1dSq4QiRVh9//DENGphXuomMjGT06NFWjkjkJdHR0QwePJiZM81LZi5YsIASJUpYOSqRFSQZFCI5EZwKm7cLlgY7Z+vGI0Qaubq6sn37drp27QrAypUruXXrlpWjEnmBoigMGTKE2bNnAxAQEEDz5s2tHJXIKpIMivzpbhCseR8ib5knj279FXRdAG+ssnZkQqSLnZ0dM2bMUPd9fHzo0aMHBw8etGJUIjc7d+4cAQEBLFq0CIABAwawevVqGTmch0kyKPKnjR/ByWUwvTKEnACtzjxopGBJa0cmRLp5e3szZcoUdf/XX3/lxRdfpH379oSHh1sxMpEbvf7662zduhWAadOmMWfOHBwdHa0clchKkgyK/OnRqiLV34Ci/lYNRYjMMHLkSCIiIpg9ezZ6vR6ADRs2UL9+fcLCwqwcncgtJk6cyOnTpwH44YcfGDp0qHUDEtlCkkGR/xyYCzEh5u0WE0BufYg8omDBggwcOJCoqChat24NwJUrV2TJMJEmUVFRfPbZZwDUrl2b/v37y63hfEKSQZG/nP4dNn9s3nYuAo6FrBuPEFnA3t6ejRs30q9fPwCOHDlCcnKylaMSOd2lS5dQFAWAbdu2WTkakZ0kGRT5R8QVWDfo8f47f0uroMizNBoNP/30E3q9nuTkZLp27ap+0Avxb1FRUYwYMQKARo0a4erqauWIRHaSZFDkH4VKwyd3oPEIeGsruMnyXSJv02q11K9fH4B169bJpNQiVdOmTcPd3Z2dO3cC8OKLL1o3IJHtJBkUeV/4BTj8PzAZQauFZmOgRB1rRyVEtti4caO6XaVKFZYsWUJSUpIVIxI5iclkYtq0aer+nDlz+PLLL60YkbAGSQZF3pacBCvegC2fwup3ISHK2hEJka0cHBzUEaEJCQn06tWLr7/+2rpBiRxj3759BAcH4+joSEJCAgMGDMDGxsbaYYlsJsmgyNsWvQwRlyA5ATzLg530gxH5zzfffMPChQvV/XHjxhEREWHFiERO8fPPPwPw2muvYWdnZ+VohLVIMijyroRouH3YvO3fExp9JANGRL6k0+no06cPO3bsAMy3BqtUqcKoUaNYvny5DCzJp7Zs2cL//vc/AHr37m3laIQ1STIo8iZDAvzYxLzt5AEdv5dEUOR7jRo14o033gAgJCSEKVOm0L17d/bt22flyER2io6OpmfPnupclO+//z5Nmza1blDCqiQZFHnTtV1w/6p5u1RTq4YiRE6h0+lYsmQJZ86cseg32L17d65fv269wES2+uGHH1i2bBkAdevWtRhAIvInSQZF3nFtDxxdbN4u6m9eXaRyV+g416phCZHTVKpUiY8++oi5c83/b9y+fVvdFnnf/v37Aejbty/79+/H3t7eyhEJa5NkUOQNQZtgcTs4utA8hYyLFzQcCl3ng05GxgmRmvfee4/OnTsDsGrVKmJjY60ckcgOhw+b+1K/9dZbaLWSBghJBkVuFx0CPzSAX1837985bt14hMhl3n33XQAuX74sgwjygTt37hASEoJWq8Xf39/a4YgcQpJBkbtd3ARhZ8zbOj303gBanXVjEiIXadmyJR07dgRg9erVVK5cmQ0bNlg3KJFltmzZApi7Cjg5OVk5GpFTSDIocreYUPPP4rXh01Ao2ci68QiRy2i1WlavXq0uW3f27Fnat2/PTz/9JFPO5CGPVhrp168fAC+99JKVIxI5iSSDIndLigOtLZRpIS2CQmSQRqNh586dfP/99+qxd999lwYNGnDp0iUrRiYyQ2xsLBUrVuTDDz8EwM/Pj48++sjKUYmcRJJBkXvEhMFPzeGbsnD6d/OxpqPA2RMKv2Dd2ITI5WxtbXn//fe5cOECer0egMDAQMqWLcvs2bOtHJ14HosWLSIoKAiAgQMHcuHCBUqUKGHlqEROIsmgyB1uHYZvy0LwEYgNA5ei5uN2LvDSJ1DhFevGJ0QeUa5cOSIjIy1ajgYPHkxUlKzrnRvdvn2bwYMHAzB9+nRmz54ty86JFCQZFDnb7q/hcw+Y3+LxsUqdoEiVx/vV3wCdbfbHJkQe5eDgwNdff622JgG4u7sTGRlpvaBEhnz11Vfq9ltvvWXFSEROJsmgyDkMCXB0EazpD3eDQFFAowNj0uMy9QZBp3lg72q1MIXIL8qWLcuYMWPU/QMHDlgxGpFeJpOJNWvWAPDll1/i7Oxs5YhETiWz8QrrSk6Ci5vB3g3WvAcxIebjtd4yryVcb5B5ypiH96FyF/CqZN14hchnPv/8c44dO8bGjRvZsGEDAQEBaGSd71zh8OHD3LlzBxcXF4YNG2btcEQOJsmgsJ7z62HFG6mfK1jK/NNGD/UHZV9MQogUfH19AZgzZw4BAQG0b9/eyhGJZ0lKSmL37t3qhOJt27aVJefEM8ltYpH9zq6B+9fM08G8/I35mOafaWFq9YPxkeBUyGrhCSEsvfPOO+r2pk2brBiJSM1vv/3G22+/zbvvvkuxYsWwt7enZcuWXLt2DYDu3btbOUKR00nLoMg+hoewcQQcXwIjr4GtA5R6Cd7eDsVrQUIU2Lmabw8LIXKM6tWr88MPP/D+++/zww8/UKxYMT799FNrh5XvJSYm0rJlS/bs2ZPinLOzM/Xr1+eTTz6hSZMmVohO5CbSMiiyz5GF5kRQ7wLbJ5iPFS5jTgTB3G9QEkEhcqSAgAB1e8yYMaxdu9Z6wQji4uIoX768RSI4YsQIFi9ezKVLl4iKimLLli2SCKZizpw5+Pn5YW9vT926dTl06NBTyxoMBiZOnEjp0qWxt7enWrVqbN682aKMn58fGo0mxWPgwIFqmaZNm6Y4379//yx7j+klLYMi64WcggPfw8lfzfuKERoOt25MQoh0KVmyJLGxseqI1H79+tGhQwcZTGIFZ86cYciQIVy/fh2ATp06sXLlSnQ6WYXpv6xYsYLhw4czd+5c6taty/Tp0wkICCAoKAhPT88U5ceMGcMvv/zCTz/9RPny5dmyZQudOnVi//79VK9eHTAP1DEajepzzpw5Q8uWLXn11VctrvXOO+8wceJEdd/R0TGL3mX6aRRZfJLbt29TokQJbt26RfHixa0dTu6w5VM4sgDcfaCoP7gVNy8NF3EJ4iPgzTXgUMA8PczEgqCYHj/37b+heE2rhS6EyLh9+/bRsGFDANzc3Ojfvz+Ojo7ExsbSqlUrWrRo8R9XEBm1dOlSPvvsMy5fvqwe+/zzzy2m/8lv0vv5XbduXWrXrq2uqmMymShRogSDBw9m1KhRKcp7e3vz6aefWrTydenSBQcHB3755ZdUX2Po0KFs2LCBS5cuqV+WmjZtir+/P9OnT8/Au8x60jIoMqbVJPMt3f2z4O6FlOc3jjTPB6iYoHRzuLwVWn4O1V43Lx8nhMiVGjRoQEBAAFu2bCEqKoopU6ao577++mv8/Pxo27YtSUlJjBo1ilKlSlkx2rwjIiKCt99+m4SEBABq1qzJZ599Rrt27awcWe6RlJTE0aNHGT16tHpMq9XSokULAgMDU31OYmJiipHYDg4O7N2796mv8csvvzB8+PAUreZLly7ll19+oUiRIrRv356xY8fmmNZBSQZFxtXsC36NzLeBY0MhMQZcvUFnB4VKm5NFrQ34NTQngVW6WjtiIUQmWLt2LUuWLCEwMBCNRoONjQ0//vgjANevX2fOnDkA/PTTT6xdu5YOHTpYM9w8YcGCBSQkJODj48OxY8coVEhmXHhSTEwM0dHR6r6dnV2KZffu3buH0WjEy8vL4riXlxcXLqTSqIG5r+y0adNo3LgxpUuXZvv27axevdritvCT1q5dS2RkJH369LE43qNHD3x9ffH29ubUqVN8/PHHBAUFsXr16gy828wnt4mR28TPFHsXzq6GyJvg7AUH54JnBajRGyrKesBCCLPo6GjmzZvHyZMnCQsLY9u2beq5+vXr07p1ax4+fEhERATdu3enadOm1gs2lzEYDJQtW5br168zf/58+vXrl6HrKIqCgpLip0kxpThuwoSjjSM2WnObUbwhnuikaBxsHHCzcwPAaDJyI+YGJpMJEyYURcGoGNWfJsWkPqp6VEWv02danTzy6PP738aPH89nn31mcezOnTsUK1aM/fv3U69ePfX4yJEj2bVrFwcPHkxxnbt37/LOO++wfv16NBoNpUuXpkWLFixYsICHDx+mKB8QEIBer2f9+vXPjPvvv/+mefPmXL58mdKlS6fx3WYdaRkUTxcdAt9VMg/4sDgeDCXqApIMCiHMXF1dGTFihLq/detW2rVrR1JSEvv372f//v3quR9//JE//vgjz09ebTAaiDPEEZ0UTawhlvIFy6PVmCfxOBp2lAv3L1C5cGWqeVQD4ErkFaYenorBZCDJmER8cjxGk5GQwyFcv34dvaueZXbLWLZy2ePEC4Xf2v2Gl5O5tWva0WksPbeU3pV680GNDwC4FX2Ll9e8nO74V7RbQcVCFQFYdmEZM47NoFOZTkxsYB4E8TD5IR3Wpq3Vd/ur2/F0zLouQufOnaNYsWLq/r9bBQEKFy6MTqcjLCzM4nhYWBhFihRJ9boeHh6sXbuWhIQEIiIi8Pb2fmr3hxs3brBt27Y0tfbVrVsXQJJBkQskxpjXAH74AHwbmPfjI6BMc6j9trWjE0LkYC1btuTevXv88ssv7N+/Hzs7O3Q6nXo7uXPnzty/fx8XFxcrR/psQfeDCI0LpZR7KUq4mFugbsfcZkXQCmINsSQZk0gyJhFriCXeEE+sIZbopGjiDHHEJMVYXOtAjwM42ToBsObSGtZdWcewmsPUZDDeEM/+O/v5t9BToQA4VnMkLDkMki3PG5/4wm40GUkyJZGsPFEogwO+n7xxqNPosNHaqMksgFajxVXvilajTfHQaXQpjmUlFxcXXF2fvWa9Xq+nZs2abN++nY4dOwLmASTbt29n0KBnr3Rlb29PsWLFMBgMrFq1itdeey1FmYULF+Lp6Unbtm3/M94TJ04AULRo0f8smx3kNjFym/iZ5tSFlz6VW8JCiEwRGBhI/fr1AXB3d+fKlSsULFjQokx0dDSbNm2iRYsWz9U3LjopmsOhhwmJDeFh8kOSlWSMJiNGxYjRZFT345PjiTPEERoXSkxSDGs7rEWnNU/T8uHOD/nrxl+MqjOKnhV6AnAk9Ah9t/RNcxwONg646F1Y1X4V7vbuACy/sJyDIQfp9EInGhdvDEBkQiR7gvdgq7NFr9Vjb2OPrdaW97u8z9H9Rxk3bRxd3+iKVqtFy+Oky9fVF1udLQBRiVE8TH6Io60jrnpzcmQ0GXmQ+ACtRovm0X//zHWnRWvefnQMjVpOp9VleRL3vNL7+b1ixQp69+7NvHnzqFOnDtOnT+e3337jwoULeHl50atXL4oVK8bkyZMBOHjwIMHBwfj7+xMcHMxnn33GtWvXOHbsGO7u7up1TSYTJUuWpHv37nz11VcWr3nlyhWWLVvGyy+/TKFChTh16hTDhg2jePHi7Nq1K1PrI6OkZVCYJUTDjf2w+2tAgVcXmaeNGZiyD4UQQmRUvXr1GDx4MLNmzSIyMpIaNWowZcoUXn75ZVxcXLhy5Qr169cnPDwcME/j8eV3X+Lo7kiSMYkEYwLh8eFcj7pOjwo91H5oS84tYf2V9bxT9R1a+rYE4HzEeYbuGJruGOOT43HRm1ssS7mXonJsZbWfHICXkxdvVHgDNzs37HR26HV6HG0ccbJ1wsnWCTc7NxxtHSloVxBnvbPa7+5Jr5d/ndfLv25xzN3enfalLW+dx8bGcnT/UQBea/UalTwqPTN2Nzs3i1gBdFodhR0Kp70C8rBu3bpx9+5dxo0bR2hoKP7+/mzevFkdVHLz5k202scJcEJCAmPGjOHq1as4Ozvz8ssvs2TJEotEEGDbtm3cvHkz1f6cer2ebdu2MX36dOLi4ihRogRdunTJUVMCScsg0jLI9b2w6F/N2g4FYdgZ0DtZJyYhRK5nUkwkJCfgaPt4+oxbMbe4/OAyP4z5gZVLVlqUdy3hSvSt6H9fBgC3em4UDiiM1kGLbQFbtHote7rtUVvaPg/8nN8u/kb/av0Z6G+eE+58xHnG7x+Pt7M37nbu6DQ6dFqd+ecT2w42DjjaOuLh4EEhh0JU9aiKrdY2ayolHWbNmsUHH5j7/bm5uXH//n2LREXI53dmkZbB/O73fnBm1eP9AiXhhVbQdJQkgkIICw8SHrAneA/Xo65zO/Y2BqOBJFOS+vNh8kOSjEkkGhOJM8TxIOEBCgoHexxUE8KFZxay8uJK+g/sT5UXqrBlyxb27dsHoCaCGhsNxfoV4/6u+8QHxQMQFRhFVGAUAC7FXej9Q2+LxKhtqbbULlqbKoWrqMcqFKrAb+1/y5a6ySxxcXGMGzeO48ePs2PHDgB0Oh0ffPCBJIIiy0gymF88fAAhJ8FkhOK1zQNDAB7ceFzmzbVQ+iWrhCeEsA5FUYhOiuZu/F2KuxTH3sY8we7u27s5EHKAioUq0q6UeWLj2zG3+XTvp+l+jXsP7+Fj6wOAh4MHFQpWwLugNwPHDmTs2LH8sfkPZsydga2tLa07t8a/pj/u7u442zpz5+odvhzzJTt37CQxMRGAmNsxzG4/mxoLalChQgUuX75MvXr1qFG6RibVivVMmDCBadOmqftt27ZVpzURIqvIbWLyQTPzrUMwv+XjfXdfaD7OPAl0yEnzKGGHguBV0XoxCiEyjaIoBMcGc+/hPcLiw7gbf5eQuBDiDHFEJERw/+F9kkxJJCQnEBYfxsNk83xpy9sup1Jhc5+0RWcW8e3Rb2lXqh2TG5k704fGhfLJ3k/wcfGhuEtxXGxd0Ov02Ght0Ov0ONg4qH3o7HR2eDp64mTrhL3OPlOSGUVReO+99/jpp59SPR8QEMCmTZtybeJ07tw5/P39MRgMODs706lTJ6ZMmZJjRpzmRHn+8zubSMtgbpYUD5E3zD+jgyEmFKq/Afp/+ufs+NKcCF7d8fg5BUuBY2HwNi+wTdFq2R+3ECLNkoxJRCZGmh8JkdyMuUnDYg0p4mSeF239lfX8eOpH6nnX45O6nwCQbEqmzeo26XodV70rsYZYdb9u0br0q9yPCoUqqMeKOBVhQcCCTHhXGaPRaPjxxx/p1asXI0eOVDv737p1C4AtW7YwZswYJk2alG0JYUhICHPmzCE4OBiNRkNERAQFChSgcOHCVKpUiXr16lGgQAECAwNxcHBg9+7d/P7779SpU4d69epRtmxZfHx8cHR0pG/fvhgMBtq2bcuGDRuyJX4hQJLB3MlkhG3jzesC/1ulTo+TwYcPLBPBVxdDpY7ZEqIQ+YXBaCDWEEucIU59xCTFcPfhXRoXb6xOtHsk9Aibrm2iXMFyvFbOPEeZoii8v+19kk3JJCvJJJvMU50kGBOISYohJimG+OT4FK85p/kcNRlMNCZyPfo6fq5+6nkbrQ0ONg4UsCuAp6Mnno6eFHEqgqveFSdbJ4o6FcXexh69To+r3hVfV1/19vAjFQpVsEgEc5KGDRtaTGIdGxtL6dKlCQ8P58svv+SPP/5gy5YteHt7Z2kcGzduTNOccqm5ePEiv/zyS6rnnlzvWYjskGeSwTlz5vD1118TGhpKtWrVmDVrFnXq1LF2WJkjORFs/plNPSEKFrSB8LOPz7v5gIM7FPADmyeW+6nV75+WwEJQ6iVw9sjOqIXI0UyKiVhDLDFJMRR2KIydzvz/2PmI8+y4tQNvZ286lukImBOuwdsHk2RKUgdIJBmTiDfEE/4w/Kmv8WPLH9Vk8GrUVX67+BvNSjRTk0GNRsOBkAMWkwanRqvR4m7njrudOwXtC6rTngA0Lt6YhQELLVZ30Gg0HOp5KEP1khs5OzuzevVqGjZsCMCZM2d45ZVXWL16NV5eXqmuRpFRMTExbNu2jT179jBrVipfyJ/Bz8+P69evP7NM7969qVTp2dPHCJHZ8kQyuGLFCoYPH87cuXOpW7cu06dPJyAggKCgIDw9s275m0ynKBC0Ea7ugtgwiA1//HP0LdBowN4NHP+ZoLXUS9BpHrh4pX49zwrmhxB5kKIoJJmSMJqM6kjVJGMSKy+uxGA0kGhMJCQuhMjESBKMCSQmJ6pLg0UnRRObFIuCucv0z21+prqnuevEsfBj/HDyB9qWaqsmg1q0BIYEPjMeBxsHHG0ccdY742TrhIeDhzrpL0ClwpUYUG0AJd1KWjzv8wafo9VosdHaYKOxUfvfuepdcdY7427njove5amT/z5q+cvvGjRoQHR0NO3atWP37t0cPXoUX19fXFxcOHv2rLp+7alTp9i4cSNHjhzh0KFD3L59G0VRaNSoEfPmzaNChaf/zbx//z7+/v7qbelHdu3aRePGjS2OKYpCfHw827Ztw8bGhlq1aqlz2YF5kuLw8HAiIiK4fv06CQkJVK5cOdVlzoTIanliAEndunWpXbs2s2fPBsz/k5UoUYLBgwczatSo/3x+juiAeusQLOsGD++nfv6jS+D8zx/8oM3gVhy8KpkTRCGygaIoGEwGDCaD+bbmo8c/qzgYTAa1xaxsgbJqgnbpwSUu3L9ACZcS+Hv6A+Y1TRedXYRJMWE0WS5o/2iB+0cL3htMBhKSE0g0JjKo+iBeKPACAEvPL+WrQ19ZDHCIN8RTd1nddL0ve509M5rNoL63eVWM/cH72XpzK3WL1qW1X2v1vW+4ugFbnS32OvPt1UcrRBR1KoqbnVuqEwuL7JeUlMSbb77Jtm3buH//8d/TWbNmce3aNWbOnElycvIzrgDt27dn9OjR1KhRQ21VVBSFt99+mwULHveZrFixItu3b3/qurYi6+WIz+88INf/9UpKSuLo0aOMHj1aPabVamnRogWBgc/+Jp8dwoMPY0iMAmMyitEAyfEQdw8lJhTKtQF387dVEiLAEA02NriUao576ebg7InBsTAhOg0aUwIlHl20XGvC4sJ4GH1DbdkALLYfbT55zEXvorYgGE1GrkVdA8wz7D9qdQiNCyU6yTzX17O+Jzx5XSdbJ3XNToAL9y9gUky84P6CukRSaFwoEQ8jUsb5r9d58pyjjSNlCpRR989GnMVgNFgkGqFxoYTEhTxetF1RMGFOKhRFMf9EsUgu7G3saVisoXrd3bd38yDhAS8WfVFd7P1q5FUOhx5WF4J/lKg8ef0nr2lSTNhqbXmv2nvqdVdeXMm1qGu0LdlWHaEZdD+I3y/+rr5PRVFQ//unHh5tP/nzs/qfqZPgLr+wnGPhx2hXqp26jNX1qOvMPD4Tk2LCYDKoS249qhc12cKEyWRS39fcFnMp5GBe7mvmsZmsvLiSNyq8ob6PG9E36Li2o1o+rX5v/zvlCpYDYNvNbXx/4nteLfuqmgwmGZP4/sT3ab7eI93Kd1OTQQcbBwDiDHHqeb1OT4BfAHqtHr1Oj5udG95O3tjZ2GGvs8fJ1gkXvQuudq646s2PRytYPFK/WH3qF6tvcUyj0aRYGULkTHq9nhUrVqAoCt9++y0jRowAYPDgwWqZli1bUr9+fZydnSlWrBjz5s2zWBZs/fr1rF+/HjBP9uzt7c358+fV89u2bSMuLo7GjRunWIlCiNwo1yeD9+7dw2g0WjS/A3h5eXHhwoVUn5OYmKjOVwXmPiBZZdCWfpzXPeVk8FrL/RLmzs79ytVjWE3zh3Fo9C3arnkZJ1snDvQ4oBYdv388++7sS1csncp0YmKDiYC5ZabTH50AOPLGEbW/1IxjM9hwNX2j2JoWb8qs5o/7zvT8sydJpiS2dt2qdnL/+dzPLDm3JF3XrVq4KkvbLlX3h/w9hLD4MFa0W0HFQuZpcDZc3cCMYzPSdd3izsXZ1GWTuj/r+Cwu3L/A3BZz1WTwePhxJh2clK7rOts6WySD225sY/+d/VQoWEFNBoNjg1ketDxd1wUYV2+cun3i7gk2XdtEpUKVHq9pmhjJ1htb031dg8mgbicYE4hMjLQYsKBBY7ng/b88Wrxep9Fhq7M1Tyui1VuM5CzpWpL63vUtbo/a6ezoWrarupj9v38+udj9o8EQep2eUm6Pb6G18G1Bc5/mONs6q8dstDZ80+SbdNeDyHs0Gg0fffQRTZo0oWfPnuh0Ojw9PenVqxf9+vWz+B3t3r07kZGRfP/990ydOpWoqCj1XFRUlMX+vHnzaN68eba+FyGyWq5PBjNi8uTJTJgwIVteS2/jiL0pDjD/4dHAP7d2NeafWhv4p1VO80+ZJ5dB0mg0ONs642jjaHHdRwufq+XQWDznSY/OPWpJeVSmoL3l4vBgTmgK2ReyKJfatZ7cdrVztSjj6eiJwWSwKOuqd6WoU9EUz33a6wB4OFoOePF29kav01vUj6veFR8XH3WBdZ1GZ158/Z9k4tGi61qNeTF2nUaHh4Pldf09/CnkUAh3O3f1WDGXYrTwaaE+R73mP4u6/zth0Wq0akL9SIBfAOULlqe0e2n1mJ+bH+9Vfe/xwvBozL8KWC4U/6g+Hm3rNI+/UbQr1Y5KhSpR06umeqy4S3E+qfsJWrTY6mwtkqp/J1hPHn/yPfet1JcuL3SxOObt7M3WrlvVa9jqbLHR2GCrtU3zIvatS7amdcnWFsfsbewZX2/8fz73WZ7sjyfE09SuXZuLFy/+Zzl3d3c++eQTPvnkE8LDwwkMDOSLL77g8OHDAJQpU4alS5fmnYGJQjwh1/cZTEpKwtHRkd9//52OHTuqx3v37k1kZCTr1q1L8Zx/twwGBwdTsWJF6XMghBBC5CLSZzBz5PqFDvV6PTVr1mT79u3qMZPJxPbt26lXr16qz7Gzs8PV1VV9uLi4pFpOCCGEECKvyxO3iYcPH07v3r2pVasWderUYfr06cTFxdG3b19rhyaEEEIIkaPliWSwW7du3L17l3HjxhEaGoq/vz+bN29OMahECCGEEEJYyhPJIMCgQYMYNGiQtcMQQgghhMhVcn2fQSGEEEIIkXGSDAohhBBC5GOSDAohhBBC5GOSDAohhBBC5GOSDAohhBBC5GOSDAohhBBC5GOSDAohhBBC5GOSDAohhBBC5GOSDAohhBBC5GOSDAohhBBC5GN5Zjm652EymQAICQmxciRCCCGESKtHn9uPPsdFxkgyCISFhQFQp04dK0cihBBCiPQKCwvDx8fH2mHkWhpFURRrB2FtycnJHD9+HC8vL7TazLtzHhMTQ8WKFTl37hwuLi6Zdt28ROoo7aSu0kfqK/2kztJO6irtsrKuTCYTYWFhVK9eHRsbad/KKEkGs1B0dDRubm5ERUXh6upq7XByJKmjtJO6Sh+pr/STOks7qau0k7rK+WQAiRBCCCFEPibJoBBCCCFEPibJYBays7Nj/Pjx2NnZWTuUHEvqKO2krtJH6iv9pM7STuoq7aSucj7pMyiEEEIIkY9Jy6AQQgghRD4myaAQQgghRD4myaAQQgghRD6W75LB3bt30759e7y9vdFoNKxdu9bifFhYGH369MHb2xtHR0dat27NpUuXLMpcuXKFTp064eHhgaurK6+99pq6iskjx44do2XLlri7u1OoUCHeffddYmNj/zO+U6dO0ahRI+zt7SlRogRTp061OH/27Fm6dOmCn58fGo2G6dOnZ6geniW319Hq1aupVasW7u7uODk54e/vz5IlSzJWGWmQ2+tr0aJFaDQai4e9vX3GKuM/5Pa6atq0aYq60mg0tG3bNmMVkga5vc4MBgMTJ06kdOnS2NvbU61aNTZv3pyxyvgPObmuEhIS6NOnD1WqVMHGxoaOHTumKBMSEkKPHj0oW7YsWq2WoUOHZqQa0mzy5MnUrl0bFxcXPD096dixI0FBQSniHjhwIIUKFcLZ2ZkuXbqkqI+bN2/Stm1bHB0d8fT0ZMSIESQnJ1uU2blzJzVq1MDOzo4yZcqwaNGi/4xPURTGjRtH0aJFcXBwoEWLFin+vb744gvq16+Po6Mj7u7uGaoHkQ+Twbi4OKpVq8acOXNSnFMUhY4dO3L16lXWrVvH8ePH8fX1pUWLFsTFxanPb9WqFRqNhr///pt9+/aRlJRE+/bt1bUR79y5Q4sWLShTpgwHDx5k8+bNnD17lj59+jwztujoaFq1aoWvry9Hjx7l66+/5rPPPuPHH39Uy8THx1OqVCm++uorihQpknkV84TcXkcFCxbk008/JTAwkFOnTtG3b1/69u3Lli1bMq+SnpDb6wvA1dWVkJAQ9XHjxo3MqZx/ye11tXr1aot6OnPmDDqdjldffTXzKulfcnudjRkzhnnz5jFr1izOnTtH//796dSpE8ePH8+8SvpHTq4ro9GIg4MDH3zwAS1atEi1TGJiIh4eHowZM4Zq1ao9X2Wkwa5duxg4cCAHDhxg69atGAwGWrVqpdYHwLBhw1i/fj0rV65k165d3Llzh86dO1u8r7Zt25KUlMT+/ftZvHgxixYtYty4cWqZa9eu0bZtW1566SVOnDjB0KFDefvtt//zb/LUqVOZOXMmc+fO5eDBgzg5OREQEEBCQoJaJikpiVdffZX3338/E2smH1LyMUBZs2aNuh8UFKQAypkzZ9RjRqNR8fDwUH766SdFURRly5YtilarVaKiotQykZGRikajUbZu3aooiqLMmzdP8fT0VIxGo1rm1KlTCqBcunTpqfF8//33SoECBZTExET12Mcff6yUK1cu1fK+vr7Kd999l673nF65vY4eqV69ujJmzJi0vennkBvra+HChYqbm1uG33NG5ca6+rfvvvtOcXFxUWJjY9P+xp9DbqyzokWLKrNnz7Z4XufOnZWePXum892nT06rqyf17t1b6dChwzPLNGnSRBkyZEiarpdZwsPDFUDZtWuXoijm925ra6usXLlSLXP+/HkFUAIDAxVFUZSNGzcqWq1WCQ0NVcv88MMPiqurq/p7MXLkSKVSpUoWr9WtWzclICDgqbGYTCalSJEiytdff60ei4yMVOzs7JRff/01RXlr/R3LK/Jdy+CzJCYmAljcItNqtdjZ2bF37161jEajsZgvyd7eHq1Wa1FGr9dbrHPs4OAAoJZJTWBgII0bN0av16vHAgICCAoK4sGDB5nwDp9fbqsjRVHYvn07QUFBNG7cOCNv+bnklvqKjY3F19eXEiVK0KFDB86ePfs8bztDcktdPWn+/Pm8/vrrODk5pfftZorcUGeJiYkpuh04ODg887pZwdp1lRtERUUB5rsrAEePHsVgMFi0ZJYvXx4fHx8CAwMB8+9AlSpV8PLyUssEBAQQHR2t/h0JDAxM0RoaEBCgXiM1165dIzQ01OJ5bm5u1K1b95nPExkjyeATHv2Sjx49mgcPHpCUlMSUKVO4ffs2ISEhALz44os4OTnx8ccfEx8fT1xcHB999BFGo1Et06xZM0JDQ/n6669JSkriwYMHjBo1CkAtk5rQ0FCL/6EAdT80NDQr3nK65ZY6ioqKwtnZGb1eT9u2bZk1axYtW7bM1LpIi9xQX+XKlWPBggWsW7eOX375BZPJRP369bl9+3am18ez5Ia6etKhQ4c4c+YMb7/9dqa8/4zIDXUWEBDAtGnTuHTpEiaTia1bt6q327OTtesqpzOZTAwdOpQGDRpQuXJlwPxvqNfrU/TF8/LyUv990/I78LQy0dHRPHz4MNV4Hj03tefllM/DvESSwSfY2tqyevVqLl68SMGCBXF0dGTHjh20adNG/Rbo4eHBypUrWb9+Pc7Ozri5uREZGUmNGjXUMpUqVWLx4sV8++23ODo6UqRIEUqWLImXl5dFGWdnZ5ydnWnTpo3V3nN65ZY6cnFx4cSJExw+fJgvvviC4cOHs3Pnzkyti7TIDfVVr149evXqhb+/P02aNGH16tV4eHgwb968zK+QZ8gNdfWk+fPnU6VKFerUqZM5FZABuaHOZsyYwQsvvED58uXR6/UMGjSIvn37WrSsZYfcUFfWNHDgQM6cOcPy5cuz/bWXLl2q1pezszN79uzJ9hjyOxtrB5DT1KxZkxMnThAVFUVSUhIeHh7UrVuXWrVqqWVatWrFlStXuHfvHjY2Nri7u1OkSBFKlSqllunRowc9evQgLCwMJycnNBoN06ZNU8ts3LgRg8EAPL7FUKRIkRSjtB7tZ9VgkYzIDXWk1WopU6YMAP7+/pw/f57JkyfTtGnTzK+Q/5Ab6utJtra2VK9encuXL2deJaRRbqmruLg4li9fzsSJEzO/EtIpp9eZh4cHa9euJSEhgYiICLy9vRk1apTFa2cXa9ZVTjZo0CA2bNjA7t27KV68uHq8SJEiJCUlERkZadE6GBYWpv77FilShEOHDllc79+/A0/7PXF1dcXBwYFXXnmFunXrqueKFSumtrKGhYVRtGhRi+f5+/s//5sWlqzdadGa+FcH49RcvHhR0Wq1ypYtW55aZvv27YpGo1EuXLjw1DLz589XHB0dlQcPHjy1zKPO2ElJSeqx0aNH56gBJKnJyXX0SN++fZUmTZo8s0xmyAv1lZycrJQrV04ZNmzYM9/H88rNdbVw4ULFzs5OuXfv3jPjz2y5uc4eSUpKUkqXLq2MHj36me/jeeW0unpSThlAYjKZlIEDByre3t7KxYsXU5x/NIDk999/V49duHAh1QEkYWFhapl58+Yprq6uSkJCgqIo5gEklStXtrh29+7d0zSA5JtvvlGPRUVFyQCSLJLvksGYmBjl+PHjyvHjxxVAmTZtmnL8+HHlxo0biqIoym+//abs2LFDuXLlirJ27VrF19dX6dy5s8U1FixYoAQGBiqXL19WlixZohQsWFAZPny4RZlZs2YpR48eVYKCgpTZs2crDg4OyowZM54ZW2RkpOLl5aW8+eabypkzZ5Tly5crjo6Oyrx589QyiYmJavxFixZVPvroI+X48eNpHsWWFrm9jr788kvlr7/+Uq5cuaKcO3dO+eabbxQbGxt1xGBmy+31NWHCBGXLli3KlStXlKNHjyqvv/66Ym9vr5w9ezaTauix3F5XjzRs2FDp1q3bc9ZG2uT2Ojtw4ICyatUq5cqVK8ru3buVZs2aKSVLlkxz4pQeObmuFEVRzp49qxw/flxp37690rRpUzXWJz06VrNmTaVHjx7K8ePHs+T/RUVRlPfff19xc3NTdu7cqYSEhKiP+Ph4tUz//v0VHx8f5e+//1aOHDmi1KtXT6lXr556Pjk5WalcubLSqlUr5cSJE8rmzZsVDw8Pi2T/6tWriqOjozJixAjl/Pnzypw5cxSdTqds3rz5mfF99dVXiru7u7Ju3Trl1KlTSocOHZSSJUsqDx8+VMvcuHFDOX78uDJhwgTF2dlZrb+YmJhMrKm8L98lgzt27FCAFI/evXsriqIoM2bMUIoXL67Y2toqPj4+ypgxYyymTVAU89QJXl5eiq2trfLCCy8o3377rWIymSzKvPnmm0rBggUVvV6vVK1aVfn555/TFN/JkyeVhg0bKnZ2dkqxYsWUr776yuL8tWvXUo0/M1u9cnsdffrpp0qZMmUUe3t7pUCBAkq9evWU5cuXZ7xC/kNur6+hQ4cqPj4+il6vV7y8vJSXX35ZOXbsWMYr5Blye10pyuOWkb/++itjlZBOub3Odu7cqVSoUEGxs7NTChUqpLz55ptKcHBwxivkGXJ6Xfn6+qYa35NSO+/r65vhOnmW1F4LUBYuXKiWefjwoTJgwAClQIECiqOjo9KpUyclJCTE4jrXr19X2rRpozg4OCiFCxdWPvzwQ8VgMFiU2bFjh+Lv76/o9XqlVKlSFq/xNCaTSRk7dqzi5eWl2NnZKc2bN1eCgoIsyvTu3TvV97Bjx46MVku+pFEURXme28xCCCGEECL3ktHEQgghhBD5mCSDQgghhBD5mCSDQgghhBD5mCSDQgghhBD5mCSDQgghhBD5mCSDQgghhBD5mCSDQgghhBD5mCSDQgghhBD5mCSDQoh8zc/Pj+nTp1s7DCGEsBpJBoUQ2aJPnz5oNBo0Gg22trZ4eXnRsmVLFixYgMlkSvN1Fi1ahLu7e7pf/2nPO3z4MO+++266ryeEEHmFJINCiGzTunVrQkJCuH79Ops2beKll15iyJAhtGvXjuTkZKvE5OHhgaOjo1VeWwghcgJJBoUQ2cbOzo4iRYpQrFgxatSowSeffMK6devYtGkTixYtAmDatGlUqVIFJycnSpQowYABA4iNjQVg586d9O3bl6ioKLWV8bPPPgMgMTGRjz76iGLFiuHk5ETdunXZuXPnfz7v37eJNRoN8+bNo127djg6OlKhQgUCAwO5fPkyTZs2xcnJifr163PlyhWL97Zu3Tpq1KiBvb09pUqVYsKECVZLcIUQIj0kGRRCWFWzZs2oVq0aq1evBkCr1TJz5kzOnj3L4sWL+fvvvxk5ciQA9evXZ/r06bi6uhISEkJISAgfffQRAIMGDSIwMJDly5dz6tQpXn31VVq3bs2lS5ee+bzUfP755/Tq1YsTJ05Qvnx5evTowXvvvcfo0aM5cuQIiqIwaNAgtfyePXvo1asXQ4YM4dy5c8ybN49FixbxxRdfZGHNCSFEJlGEECIb9O7dW+nQoUOq57p166ZUqFAh1XMrV65UChUqpO4vXLhQcXNzsyhz48YNRafTKcHBwRbHmzdvrowePfqpz1MURfH19VW+++47dR9QxowZo+4HBgYqgDJ//nz12K+//qrY29tbvM6XX35pcd0lS5YoRYsWTfU9CSFETmJj7WRUCCEURUGj0QCwbds2Jk+ezIULF4iOjiY5OZmEhATi4+Of2rfv9OnTGI1GypYta3E8MTGRQoUKpTueqlWrqtteXl4AVKlSxeJYQkIC0dHRuLq6cvLkSfbt22fREmg0Gv8zbiGEyAkkGRRCWN358+cpWbIk169fp127drz//vt88cUXFCxYkL179/LWW2+RlJT01KQqNjYWnU7H0aNH0el0FuecnZ3THY+tra26/ShJTe3Yo1HQsbGxTJgwgc6dO6e4lr29fbpfXwghspMkg0IIq/r77785ffo0w4YN4+jRo5hMJr799lu0WnOX5t9++82ivF6vx2g0WhyrXr06RqOR8PBwGjVqlOrrpPa8zFKjRg2CgoIoU6ZMllxfCCGykiSDQohsk5iYSGhoKEajkbCwMDZv3szkyZNp164dvXr14syZMxgMBmbNmkX79u3Zt28fc+fOtbiGn58fsbGxbN++nWrVquHo6EjZsmXp2bMnvXr14ttvv6V69ercvXuX7du3U7VqVdq2bZvq8zLr9u24ceNo164dPj4+dO3aFa1Wy8mTJzlz5gyTJk3KlNcQQoisIqOJhRDZZvPmzRQtWhQ/Pz9at27Njh07mDlzJuvWrUOn01GtWjWmTZvGlClTqFy5MkuXLmXy5MkW16hfvz79+/enW7dueHh4MHXqVAAWLlxIr169+PDDDylXrhwdO3bk8OHD+Pj4PPN5mSEgIIANGzbw119/Ubt2bV588UW+++47fH19M+01hBAiq2gURVGsHYQQQgghhLAOaRkUQgghhMjHJBkUQgghhMjHJBkUQgghhMjHJBkUQgghhMjHJBkUQgghhMjHJBkUQgghhMjHJBkUQgghhMjHJBkUQgghhMjHJBkUQgghhMjHJBkUQgghhMjHJBkUQgghhMjHJBkUQgghhMjH/g+ObD+CFgG6LQAAAABJRU5ErkJggg==", "text/plain": [ "
" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "from cycler import cycler\n", "plt.style.use('default')\n", "\n", "fig, ax = plt.subplots()\n", "\n", "ax.set_prop_cycle(cycler('color', ['tab:blue', 'tab:orange', 'tab:green']) + cycler('linestyle', ['-', '--', '-.']))\n", "\n", "ax.plot(timesteps['Datetime'], timesteps[['NA', 'NB', 'NC']].values)\n", "ax.legend(labels = ['$N_A$', '$N_B$', '$N_C$'], loc = 'upper left')\n", "ax.set_ylabel('Defect state percentages [%]')\n", "ax.set_xlabel('Datetime')\n", "\n", "ax2 = ax.twinx()\n", "ax2.plot(timesteps['Datetime'], timesteps['Pmp_norm'], c = 'black', label = 'Normalized STC $P_{MP}$')\n", "ax2.legend(loc = 'upper right')\n", "ax2.set_ylabel('Normalized STC $P_{MP}$')\n", "\n", "ax.set_title('Outdoor LETID \\n'f'{location.name}')\n", "\n", "plt.show()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The example data provided for Golden, CO, shows how $N_A$ increases in cold weather, and power temporarily recovers, due to temporary recovery of LETID (B->A)." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "##### The function `calc_letid_outdoors` wraps all of the steps above into a single function:" ] }, { "cell_type": "code", "execution_count": 23, "metadata": {}, "outputs": [ { "data": { "text/html": [ "
\n", "\n", "\n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", "
TemperatureInjectionNANBNCtauJscVocIscFFPmpPmp_norm
time
1999-01-01 00:30:00-07:000.0NaN100.0000000.000000e+000.000000e+00115.00000041.5909970.66632710.1066120.8409875.6634671.000000
1999-01-01 01:30:00-07:000.0NaN100.0000001.702422e-150.000000e+00115.00000041.5909970.66632710.1066120.8409875.6634671.000000
1999-01-01 02:30:00-07:000.0NaN100.0000003.404843e-155.403329e-36115.00000041.5909970.66632710.1066120.8409875.6634671.000000
1999-01-01 03:30:00-07:000.0NaN100.0000005.107265e-151.620999e-35115.00000041.5909970.66632710.1066120.8409875.6634671.000000
1999-01-01 04:30:00-07:000.0NaN100.0000006.809686e-153.241997e-35115.00000041.5909970.66632710.1066120.8409875.6634671.000000
.......................................
1999-12-31 19:30:00-07:000.0NaN27.8336546.778463e+014.381716e+0066.11214241.3338510.65497410.0441260.8389665.5192570.974537
1999-12-31 20:30:00-07:000.0NaN27.8336546.778463e+014.381716e+0066.11214241.3338510.65497410.0441260.8389665.5192570.974537
1999-12-31 21:30:00-07:000.0NaN27.8336546.778463e+014.381716e+0066.11214241.3338510.65497410.0441260.8389665.5192570.974537
1999-12-31 22:30:00-07:000.0NaN27.8336546.778463e+014.381716e+0066.11214241.3338510.65497410.0441260.8389665.5192570.974537
1999-12-31 23:30:00-07:000.0NaN27.8336546.778463e+014.381716e+0066.11214241.3338510.65497410.0441260.8389665.5192570.974537
\n", "

8760 rows × 12 columns

\n", "
" ], "text/plain": [ " Temperature Injection NA NB \\\n", "time \n", "1999-01-01 00:30:00-07:00 0.0 NaN 100.000000 0.000000e+00 \n", "1999-01-01 01:30:00-07:00 0.0 NaN 100.000000 1.702422e-15 \n", "1999-01-01 02:30:00-07:00 0.0 NaN 100.000000 3.404843e-15 \n", "1999-01-01 03:30:00-07:00 0.0 NaN 100.000000 5.107265e-15 \n", "1999-01-01 04:30:00-07:00 0.0 NaN 100.000000 6.809686e-15 \n", "... ... ... ... ... \n", "1999-12-31 19:30:00-07:00 0.0 NaN 27.833654 6.778463e+01 \n", "1999-12-31 20:30:00-07:00 0.0 NaN 27.833654 6.778463e+01 \n", "1999-12-31 21:30:00-07:00 0.0 NaN 27.833654 6.778463e+01 \n", "1999-12-31 22:30:00-07:00 0.0 NaN 27.833654 6.778463e+01 \n", "1999-12-31 23:30:00-07:00 0.0 NaN 27.833654 6.778463e+01 \n", "\n", " NC tau Jsc Voc \\\n", "time \n", "1999-01-01 00:30:00-07:00 0.000000e+00 115.000000 41.590997 0.666327 \n", "1999-01-01 01:30:00-07:00 0.000000e+00 115.000000 41.590997 0.666327 \n", "1999-01-01 02:30:00-07:00 5.403329e-36 115.000000 41.590997 0.666327 \n", "1999-01-01 03:30:00-07:00 1.620999e-35 115.000000 41.590997 0.666327 \n", "1999-01-01 04:30:00-07:00 3.241997e-35 115.000000 41.590997 0.666327 \n", "... ... ... ... ... \n", "1999-12-31 19:30:00-07:00 4.381716e+00 66.112142 41.333851 0.654974 \n", "1999-12-31 20:30:00-07:00 4.381716e+00 66.112142 41.333851 0.654974 \n", "1999-12-31 21:30:00-07:00 4.381716e+00 66.112142 41.333851 0.654974 \n", "1999-12-31 22:30:00-07:00 4.381716e+00 66.112142 41.333851 0.654974 \n", "1999-12-31 23:30:00-07:00 4.381716e+00 66.112142 41.333851 0.654974 \n", "\n", " Isc FF Pmp Pmp_norm \n", "time \n", "1999-01-01 00:30:00-07:00 10.106612 0.840987 5.663467 1.000000 \n", "1999-01-01 01:30:00-07:00 10.106612 0.840987 5.663467 1.000000 \n", "1999-01-01 02:30:00-07:00 10.106612 0.840987 5.663467 1.000000 \n", "1999-01-01 03:30:00-07:00 10.106612 0.840987 5.663467 1.000000 \n", "1999-01-01 04:30:00-07:00 10.106612 0.840987 5.663467 1.000000 \n", "... ... ... ... ... \n", "1999-12-31 19:30:00-07:00 10.044126 0.838966 5.519257 0.974537 \n", "1999-12-31 20:30:00-07:00 10.044126 0.838966 5.519257 0.974537 \n", "1999-12-31 21:30:00-07:00 10.044126 0.838966 5.519257 0.974537 \n", "1999-12-31 22:30:00-07:00 10.044126 0.838966 5.519257 0.974537 \n", "1999-12-31 23:30:00-07:00 10.044126 0.838966 5.519257 0.974537 \n", "\n", "[8760 rows x 12 columns]" ] }, "execution_count": 23, "metadata": {}, "output_type": "execute_result" } ], "source": [ "nA_0 = 100\n", "nB_0 = 0\n", "nC_0 = 0\n", "mechanism_params = 'repins'\n", "\n", "letid.calc_letid_outdoors(tau_0, tau_deg,wafer_thickness, s_rear, nA_0, nB_0, nC_0, weather, meta, mechanism_params, generation_df, module_parameters = cec_module)" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [] } ], "metadata": { "kernelspec": { "display_name": "Python 3 (ipykernel)", "language": "python", "name": "python3" }, "language_info": { "codemirror_mode": { "name": "ipython", "version": 3 }, "file_extension": ".py", "mimetype": "text/x-python", "name": "python", "nbconvert_exporter": "python", "pygments_lexer": "ipython3", "version": "3.11.7" }, "vscode": { "interpreter": { "hash": "848658e0671c41dd18b216771b1713479db7d685859cbb6c795b270024b1888c" } } }, "nbformat": 4, "nbformat_minor": 2 }